source: trunk/LMDZ.COMMON/libf/phy_common/vertical_layers_mod.F90

Last change on this file was 1621, checked in by emillour, 8 years ago

Further work on full dynamics/physics separation.

LMDZ.COMMON:

  • added phy_common/vertical_layers_mod.F90 to store information on vertical grid. This is where routines in the physics should get the information.
  • The contents of vertical_layers_mod intialized via dynphy_lonlat/inigeomphy_mod.F90.

LMDZ.MARS:

  • physics now completely decoupled from dynamics; the physics package may now be compiled as a library (-libphy option of makelmdz_fcm).
  • created an "ini_tracer_mod" routine in module "tracer_mod" for a cleaner initialization of the later.
  • removed some purely dynamics-related outputs (etot0, zoom parameters, etc.) from diagfi.nc and stats.nc outputs as these informations are not available in the physics.

LMDZ.GENERIC:

  • physics now completely decoupled from dynamics; the physics package may now be compiled as a library (-libphy option of makelmdz_fcm).
  • added nqtot to tracer_h.F90.
  • removed some purely dynamics-related outputs (etot0, zoom parameters, etc.) from diagfi.nc and stats.nc outputs as these informations are not available in the physics.

LMDZ.VENUS:

  • physics now completely decoupled from dynamics; the physics package may now be compiled as a library (-libphy option of makelmdz_fcm).
  • added infotrac_phy.F90 to store information on tracers in the physics. Initialized via iniphysiq.
  • added cpdet_phy_mod.F90 to store t2tpot etc. functions to be used in the physics. Initialized via iniphysiq. IMPORTANT: there are some hard-coded constants! These should match what is in cpdet_mod.F90 in the dynamics.
  • got rid of references to moyzon_mod module within the physics. The required variables (tmoy, plevmoy) are passed to the physics as arguments to physiq.

LMDZ.TITAN:

  • added infotrac_phy.F90 to store information on tracers in the physics. Initialized via iniphysiq.
  • added cpdet_phy_mod.F90 to store t2tpot etc. functions to be used in the physics.
  • Extra work required to completely decouple physics and dynamics: moyzon_mod should be cleaned up and information passed from dynamics to physics as as arguments. Likewise moyzon_ch and moyzon_mu should not be queried from logic_mod (which is in the dynamics).

EM

File size: 2.3 KB
RevLine 
[1621]1! $Id: $
2
3MODULE vertical_layers_mod
4
5   REAL,SAVE             :: preff  ! reference surface pressure (Pa)
6   REAL,SAVE             :: scaleheight ! atmospheric reference scale height (km)
7   REAL,SAVE,ALLOCATABLE :: ap(:) ! hybrid (pressure contribution) coordinate
8                                  ! at layer interfaces (Pa)
9   REAL,SAVE,ALLOCATABLE :: bp(:) ! hybrid (sigma contribution) coordinate
10                                  ! at layer interfaces
11   REAL,SAVE,ALLOCATABLE :: aps(:) ! hybrid (pressure contribution) coordinate
12                                   ! at mid-layer (Pa)
13   REAL,SAVE,ALLOCATABLE :: bps(:) ! hybrid (sigma contribution) coordinate
14                                   ! at mid-layer
15   REAL,SAVE,ALLOCATABLE :: presnivs(:) ! reference pressure at mid-layer (Pa),
16                                        ! based on preff, ap and bp
17   REAL,SAVE,ALLOCATABLE :: pseudoalt(:) ! pseudo-altitude of model layers (km),
18                                         ! based on preff and scaleheight
19   
20!$OMP THREADPRIVATE(preff,scaleheight,ap,bp,aps,bps,presnivs,pseudoalt)
21
22
23CONTAINS
24
25  SUBROUTINE init_vertical_layers(nlayer,preff_,scaleheight_,ap_,bp_,&
26                                 aps_,bps_,presnivs_, pseudoalt_)
27    IMPLICIT NONE
28    INTEGER,INTENT(IN) :: nlayer ! number of atmospheric layers
29    REAL,INTENT(IN)    :: preff_ ! reference surface pressure (Pa)
30    REAL,INTENT(IN)    :: scaleheight_ ! atmospheric scale height (km)
31    REAL,INTENT(IN)    :: ap_(nlayer+1) ! hybrid coordinate at interfaces
32    REAL,INTENT(IN)    :: bp_(nlayer+1) ! hybrid coordinate at interfaces
33    REAL,INTENT(IN)    :: aps_(nlayer) ! hybrid coordinate at mid-layer
34    REAL,INTENT(IN)    :: bps_(nlayer) ! hybrid coordinate at mid-layer
35    REAL,INTENT(IN)    :: presnivs_(nlayer) ! Appproximative pressure of atm. layers (Pa)
36    REAL,INTENT(IN)    :: pseudoalt_(nlayer) ! pseudo-altitude of atm. layers (km)
37 
38    ALLOCATE(ap(nlayer+1))
39    ALLOCATE(bp(nlayer+1))
40    ALLOCATE(aps(nlayer))
41    ALLOCATE(bps(nlayer))
42    ALLOCATE(presnivs(nlayer))
43    ALLOCATE(pseudoalt(nlayer))
44 
45    preff = preff_
46    scaleheight=scaleheight_
47    ap(:) = ap_(:)
48    bp(:) = bp_(:)
49    aps(:) = aps_(:)
50    bps(:) = bps_(:)
51    presnivs(:) = presnivs_(:)
52    pseudoalt(:) = pseudoalt_(:)
53
54  END SUBROUTINE init_vertical_layers
55
56END MODULE vertical_layers_mod
Note: See TracBrowser for help on using the repository browser.