source: trunk/LMDZ.MARS/libf/phymars/comsoil_h.F90 @ 2909

Last change on this file since 2909 was 2909, checked in by romain.vande, 21 months ago

Mars PEM:
New Boolean options for following orbital parameters of ob_ex_lsp.asc: var_obl, var_ex, var_lsp.
If using evol_orbit_pem=true, you can specify which parameter to follow.
True by default: Do you want to change the parameter XXX during the PEM run as prescribed in ob_ex_lsp.asc.
If false, it is set to constant (to the value of the tab_cntrl in the start)
RV

File size: 3.6 KB
Line 
1module comsoil_h
2
3implicit none
4! nsoilmx : number of subterranean layers
5!EM: old soil routine:      integer, parameter :: nsoilmx = 10
6  integer, parameter :: nsoilmx = 18
7
8  real,save,allocatable,dimension(:) :: layer      ! soil layer depths
9  real,save,allocatable,dimension(:) :: mlayer     ! soil mid-layer depths
10  real,save,allocatable,dimension(:,:) :: inertiedat ! soil thermal inertia
11  real,save :: volcapa    ! soil volumetric heat capacity
12       ! NB: volcapa is read fromn control(35) from physicq start file
13       !     in physdem (or set via tabfi, or initialized in
14       !                 soil_settings.F)
15
16!$OMP THREADPRIVATE(layer,mlayer,inertiedat,volcapa)
17
18  ! variables (FC: built in firstcall in soil.F)
19  REAL,SAVE,ALLOCATABLE :: tsoil(:,:,:)       ! sub-surface temperatures (K)
20  real,save,allocatable :: mthermdiff(:,:,:)  ! (FC) mid-layer thermal diffusivity
21  real,save,allocatable :: thermdiff(:,:,:)   ! (FC) inter-layer thermal diffusivity
22  real,save,allocatable :: coefq(:)         ! (FC) q_{k+1/2} coefficients
23  real,save,allocatable :: coefd(:,:,:)       ! (FC) d_k coefficients
24  real,save,allocatable :: alph(:,:,:)        ! (FC) alpha_k coefficients
25  real,save,allocatable :: beta(:,:,:)        ! beta_k coefficients
26  real,save :: mu
27  real,save,allocatable :: flux_geo(:)       ! Geothermal Flux (W/m^2)
28
29!$OMP THREADPRIVATE(tsoil,mthermdiff,thermdiff,coefq,coefd,alph,beta,mu,flux_geo)
30
31contains
32
33  subroutine ini_comsoil_h(ngrid,nslope)
34 
35  implicit none
36  integer,intent(in) :: ngrid ! number of atmospheric columns
37  integer,intent(in) :: nslope ! number of sub grid slopes
38    allocate(layer(nsoilmx)) !soil layer depths
39    allocate(mlayer(0:nsoilmx-1)) ! soil mid-layer depths
40    allocate(inertiedat(ngrid,nsoilmx)) ! soil thermal inertia
41 
42    allocate(tsoil(ngrid,nsoilmx,nslope)) ! soil temperatures
43
44    allocate(mthermdiff(ngrid,0:nsoilmx-1,nslope))
45    allocate(thermdiff(ngrid,nsoilmx-1,nslope))
46    allocate(coefq(0:nsoilmx-1))
47    allocate(coefd(ngrid,nsoilmx-1,nslope))
48    allocate(alph(ngrid,nsoilmx-1,nslope))
49    allocate(beta(ngrid,nsoilmx-1,nslope))
50    allocate(flux_geo(ngrid))
51 
52  end subroutine ini_comsoil_h
53
54
55  subroutine end_comsoil_h
56
57  implicit none
58
59    if (allocated(layer)) deallocate(layer)
60    if (allocated(mlayer)) deallocate(mlayer)
61    if (allocated(inertiedat)) deallocate(inertiedat)
62    if (allocated(tsoil)) deallocate(tsoil)
63    if (allocated(mthermdiff)) deallocate(mthermdiff)
64    if (allocated(thermdiff)) deallocate(thermdiff)
65    if (allocated(coefq)) deallocate(coefq)
66    if (allocated(coefd)) deallocate(coefd)
67    if (allocated(alph)) deallocate(alph)
68    if (allocated(beta)) deallocate(beta)
69    if (allocated(flux_geo)) deallocate(flux_geo)
70  end subroutine end_comsoil_h
71
72  subroutine ini_comsoil_h_slope_var(ngrid,nslope)
73 
74  implicit none
75  integer,intent(in) :: ngrid ! number of atmospheric columns
76  integer,intent(in) :: nslope ! number of sub grid slopes
77 
78    allocate(tsoil(ngrid,nsoilmx,nslope)) ! soil temperatures
79    allocate(mthermdiff(ngrid,0:nsoilmx-1,nslope))
80    allocate(thermdiff(ngrid,nsoilmx-1,nslope))
81    allocate(coefd(ngrid,nsoilmx-1,nslope))
82    allocate(alph(ngrid,nsoilmx-1,nslope))
83    allocate(beta(ngrid,nsoilmx-1,nslope))
84 
85  end subroutine ini_comsoil_h_slope_var
86
87
88  subroutine end_comsoil_h_slope_var
89
90  implicit none
91
92    if (allocated(tsoil)) deallocate(tsoil)
93    if (allocated(mthermdiff)) deallocate(mthermdiff)
94    if (allocated(thermdiff)) deallocate(thermdiff)
95    if (allocated(coefd)) deallocate(coefd)
96    if (allocated(alph)) deallocate(alph)
97    if (allocated(beta)) deallocate(beta)
98
99  end subroutine end_comsoil_h_slope_var
100
101end module comsoil_h
Note: See TracBrowser for help on using the repository browser.