source: trunk/LMDZ.COMMON/libf/evolution/soil_settings_PEM_mod.F90 @ 3428

Last change on this file since 3428 was 3206, checked in by jbclement, 10 months ago

PEM:
Cleanings of unused variables/arguments and bad type conversions.
JBC

File size: 3.9 KB
Line 
1MODULE soil_settings_PEM_mod
2
3implicit none
4
5!=======================================================================
6contains
7!=======================================================================
8
9SUBROUTINE soil_settings_PEM(ngrid,nslope,nsoil_PEM,nsoil_PCM,TI_PCM,TI_PEM)
10
11use comsoil_h_PEM, only: layer_PEM, mlayer_pem, depth_breccia, depth_bedrock, index_breccia, index_bedrock
12use iostart,       only: inquire_field_ndims, get_var, get_field, inquire_field, inquire_dimension_length
13
14implicit none
15
16!=======================================================================
17!  Author: LL, based on work by  Ehouarn Millour (07/2006)
18!
19!  Purpose: Read and/or initialise soil depths and properties
20!
21! Modifications: Aug.2010 EM: use NetCDF90 to load variables (enables using
22!                              r4 or r8 restarts independently of having compiled
23!                              the GCM in r4 or r8)
24!                June 2013 TN: Possibility to read files with a time axis
25!
26!  The various actions and variable read/initialized are:
27!  1. Read/build layer (and midlayer) depth
28!  2. Interpolate thermal inertia and temperature on the new grid, if necessary
29!=======================================================================
30
31!=======================================================================
32!  arguments
33!  ---------
34!  inputs:
35integer,                                 intent(in) :: ngrid     ! # of horizontal grid points
36integer,                                 intent(in) :: nslope    ! # of subslope wihtin the mesh
37integer,                                 intent(in) :: nsoil_PEM ! # of soil layers in the PEM
38integer,                                 intent(in) :: nsoil_PCM ! # of soil layers in the PCM
39real, dimension(ngrid,nsoil_PCM,nslope), intent(in) :: TI_PCM    ! Thermal inertia  in the PCM [SI]
40
41real, dimension(ngrid,nsoil_PEM,nslope), intent(inout) :: TI_PEM ! Thermal inertia   in the PEM [SI]
42!=======================================================================
43! local variables:
44integer :: ig, iloop, islope ! loop counters
45real    :: alpha, lay1       ! coefficients for building layers
46real    :: index_powerlaw    ! coefficient to match the power low grid with the exponential one
47!=======================================================================
48! 1. Depth coordinate
49! -------------------
50! mlayer_PEM distribution: For the grid in common with the PCM: lay1*(1+k^(2.9)*(1-exp(-k/20))
51! Then we use a power low : lay1*alpha**(k-1/2)
52lay1 = 2.e-4
53alpha = 2
54do iloop = 0,nsoil_PCM - 1
55    mlayer_PEM(iloop) = lay1*(1+iloop**2.9*(1-exp(-real(iloop)/20.)))
56enddo
57
58do iloop = 0, nsoil_PEM-1
59    if(lay1*(alpha**(iloop-0.5)) > mlayer_PEM(nsoil_PCM-1)) then
60        index_powerlaw = iloop
61        exit
62    endif
63enddo
64
65do iloop = nsoil_PCM, nsoil_PEM-1
66    mlayer_PEM(iloop) = lay1*(alpha**(index_powerlaw + (iloop-nsoil_PCM)-0.5))
67enddo
68
69! Build layer_PEM()
70do iloop=1,nsoil_PEM-1
71layer_PEM(iloop)=(mlayer_PEM(iloop)+mlayer_PEM(iloop-1))/2.
72enddo
73layer_PEM(nsoil_PEM)=2*mlayer_PEM(nsoil_PEM-1) - mlayer_PEM(nsoil_PEM-2)
74
75
76! 2. Thermal inertia (note: it is declared in comsoil_h)
77! ------------------
78do ig = 1,ngrid
79    do islope = 1,nslope
80        do iloop = 1,nsoil_PCM
81            TI_PEM(ig,iloop,islope) = TI_PCM(ig,iloop,islope)
82        enddo
83        if (nsoil_PEM > nsoil_PCM) then
84            do iloop = nsoil_PCM + 1,nsoil_PEM
85                TI_PEM(ig,iloop,islope) = TI_PCM(ig,nsoil_PCM,islope)
86            enddo
87        endif
88    enddo
89enddo
90
91! 3. Index for subsurface layering
92! ------------------
93index_breccia = 1
94do iloop = 1,nsoil_PEM - 1
95    if (depth_breccia >= layer_PEM(iloop)) then
96        index_breccia = iloop
97    else
98        exit
99    endif
100enddo
101
102index_bedrock = 1
103do iloop = 1,nsoil_PEM - 1
104    if (depth_bedrock >= layer_PEM(iloop)) then
105        index_bedrock = iloop
106    else
107        exit
108    endif
109enddo
110
111END SUBROUTINE soil_settings_PEM
112
113END MODULE soil_settings_PEM_mod
Note: See TracBrowser for help on using the repository browser.