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

Last change on this file since 3189 was 3149, checked in by jbclement, 2 years ago

PEM:

  • Simplification of the algorithm managing the stopping criteria;
  • Complete rework of the ice management in the PEM (H2O & CO2);

    Subroutines to evolve the H2O and CO2 ice are now in the same module "evol_ice_mod.F90".
    Tendencies are computed from the variation of "ice + frost" between the 2 PCM runs.
    Evolving ice in the PEM is now called 'h2o_ice' or 'co2_ice' (not anymore in 'qsurf' and free of 'water_reservoir').
    Default value 'ini_h2o_bigreservoir' (= 10 m) initializes the H2O ice of the first PEM run where there is 'watercap'. For the next PEM runs, initialization is done with the value kept in "startpem.nc". CO2 ice is taken from 'perennial_co2ice' of the PCM (paleoclimate flag must be true).
    Simplification of the condition to compute the surface ice cover needed for the stopping criteria.
    Frost ('qsurf') is not evolved by the PEM and given back to the PCM.
    New default threshold value 'inf_h2oice_threshold' (= 2 m) to decide at the end of the PEM run if the H2O ice should be 'watercap' or not for the next PCM runs. If H2O ice cannot be 'watercap', then the remaining H2O ice is transferred to the frost ('qsurf').

  • Renaming of variables/subroutines for clarity;
  • Some cleanings throughout the code;
  • Small updates in files of the deftank.

JBC

File size: 4.0 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
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_PCM ! # of soil layers in the PCM
42real, dimension(ngrid,nsoil_PCM,nslope), intent(in) :: TI_PCM    ! Thermal inertia  in the PCM [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
51real    :: index_powerlaw ! coefficient to match the power low grid with the exponential one
52!=======================================================================
53! 1. Depth coordinate
54! -------------------
55! mlayer_PEM distribution: For the grid in common with the PCM: lay1*(1+k^(2.9)*(1-exp(-k/20))
56! Then we use a power low : lay1*alpha**(k-1/2)
57lay1 = 2.e-4
58alpha = 2
59do iloop = 0,nsoil_PCM - 1
60    mlayer_PEM(iloop) = lay1*(1+iloop**2.9*(1-exp(-real(iloop)/20.)))
61enddo
62
63do iloop = 0, nsoil_PEM-1
64    if(lay1*(alpha**(iloop-0.5)) > mlayer_PEM(nsoil_PCM-1)) then
65        index_powerlaw = iloop
66        exit
67    endif
68enddo
69
70do iloop = nsoil_PCM, nsoil_PEM-1
71    mlayer_PEM(iloop) = lay1*(alpha**(index_powerlaw + (iloop-nsoil_PCM)-0.5))
72enddo
73
74! Build layer_PEM()
75do iloop=1,nsoil_PEM-1
76layer_PEM(iloop)=(mlayer_PEM(iloop)+mlayer_PEM(iloop-1))/2.
77enddo
78layer_PEM(nsoil_PEM)=2*mlayer_PEM(nsoil_PEM-1) - mlayer_PEM(nsoil_PEM-2)
79
80
81! 2. Thermal inertia (note: it is declared in comsoil_h)
82! ------------------
83do ig = 1,ngrid
84    do islope = 1,nslope
85        do iloop = 1,nsoil_PCM
86            TI_PEM(ig,iloop,islope) = TI_PCM(ig,iloop,islope)
87        enddo
88        if (nsoil_PEM > nsoil_PCM) then
89            do iloop = nsoil_PCM + 1,nsoil_PEM
90                TI_PEM(ig,iloop,islope) = TI_PCM(ig,nsoil_PCM,islope)
91            enddo
92        endif
93    enddo
94enddo
95
96! 3. Index for subsurface layering
97! ------------------
98index_breccia = 1
99do iloop = 1,nsoil_PEM - 1
100    if (depth_breccia >= layer_PEM(iloop)) then
101        index_breccia = iloop
102    else
103        exit
104    endif
105enddo
106
107index_bedrock = 1
108do iloop = 1,nsoil_PEM - 1
109    if (depth_bedrock >= layer_PEM(iloop)) then
110        index_bedrock = iloop
111    else
112        exit
113    endif
114enddo
115
116END SUBROUTINE soil_settings_PEM
117
118END MODULE soil_settings_PEM_mod
Note: See TracBrowser for help on using the repository browser.