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

Last change on this file since 3088 was 3076, checked in by jbclement, 21 months ago

PEM:
Big cleaning/improvements of the PEM:

  • Conversion of "abort_pem.F" and "soil_settings_PEM.F" into Fortran 90;
  • Transformation of every PEM subroutines into module;
  • Rewriting of many subroutines with modern Fortran syntax;
  • Correction of a bug in "pem.F90" when calling 'recomp_tend_co2_slope'. The arguments were given in disorder and emissivity was missing;
  • Update of "launch_pem.sh" in deftank.

JBC

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