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

Last change on this file since 2932 was 2919, checked in by llange, 22 months ago

PCM
Flux_geo is now correctled initialized (wasnot correctly read in the startfi) for the 3D and 1D
Minor fix in the pem (wrong name for a variable)
LL

File size: 3.6 KB
RevLine 
[1047]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
[2616]16!$OMP THREADPRIVATE(layer,mlayer,inertiedat,volcapa)
[2578]17
[1770]18  ! variables (FC: built in firstcall in soil.F)
[2900]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
[1770]22  real,save,allocatable :: coefq(:)         ! (FC) q_{k+1/2} coefficients
[2900]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
[1224]26  real,save :: mu
[2919]27  real,save,allocatable :: flux_geo(:,:)       ! Geothermal Flux (W/m^2)
[1224]28
[2887]29!$OMP THREADPRIVATE(tsoil,mthermdiff,thermdiff,coefq,coefd,alph,beta,mu,flux_geo)
[2578]30
[1047]31contains
32
[2900]33  subroutine ini_comsoil_h(ngrid,nslope)
[1047]34 
35  implicit none
36  integer,intent(in) :: ngrid ! number of atmospheric columns
[2900]37  integer,intent(in) :: nslope ! number of sub grid slopes
[1047]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
[1224]41 
[2900]42    allocate(tsoil(ngrid,nsoilmx,nslope)) ! soil temperatures
[1224]43
[2900]44    allocate(mthermdiff(ngrid,0:nsoilmx-1,nslope))
45    allocate(thermdiff(ngrid,nsoilmx-1,nslope))
[1224]46    allocate(coefq(0:nsoilmx-1))
[2900]47    allocate(coefd(ngrid,nsoilmx-1,nslope))
48    allocate(alph(ngrid,nsoilmx-1,nslope))
49    allocate(beta(ngrid,nsoilmx-1,nslope))
[2919]50    allocate(flux_geo(ngrid,nslope))
[1224]51 
[1047]52  end subroutine ini_comsoil_h
53
[1770]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)
[2887]69    if (allocated(flux_geo)) deallocate(flux_geo)
[1770]70  end subroutine end_comsoil_h
71
[2909]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
[1047]101end module comsoil_h
Note: See TracBrowser for help on using the repository browser.