[3076] | 1 | MODULE soil_settings_PEM_mod |
---|
[2794] | 2 | |
---|
[3076] | 3 | implicit none |
---|
[2794] | 4 | |
---|
[3076] | 5 | !======================================================================= |
---|
| 6 | contains |
---|
| 7 | !======================================================================= |
---|
| 8 | |
---|
[3149] | 9 | SUBROUTINE soil_settings_PEM(ngrid,nslope,nsoil_PEM,nsoil_PCM,TI_PCM,TI_PEM) |
---|
[3076] | 10 | |
---|
[3206] | 11 | use comsoil_h_PEM, only: layer_PEM, mlayer_pem, depth_breccia, depth_bedrock, index_breccia, index_bedrock |
---|
| 12 | use iostart, only: inquire_field_ndims, get_var, get_field, inquire_field, inquire_dimension_length |
---|
[2895] | 13 | |
---|
[3076] | 14 | implicit none |
---|
[2794] | 15 | |
---|
[3076] | 16 | !======================================================================= |
---|
[2855] | 17 | ! Author: LL, based on work by Ehouarn Millour (07/2006) |
---|
[2794] | 18 | ! |
---|
| 19 | ! Purpose: Read and/or initialise soil depths and properties |
---|
| 20 | ! |
---|
[3076] | 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 |
---|
[2794] | 25 | ! |
---|
| 26 | ! The various actions and variable read/initialized are: |
---|
[2855] | 27 | ! 1. Read/build layer (and midlayer) depth |
---|
[3076] | 28 | ! 2. Interpolate thermal inertia and temperature on the new grid, if necessary |
---|
| 29 | !======================================================================= |
---|
[2794] | 30 | |
---|
[3076] | 31 | !======================================================================= |
---|
[2794] | 32 | ! arguments |
---|
| 33 | ! --------- |
---|
| 34 | ! inputs: |
---|
[3076] | 35 | integer, intent(in) :: ngrid ! # of horizontal grid points |
---|
| 36 | integer, intent(in) :: nslope ! # of subslope wihtin the mesh |
---|
| 37 | integer, intent(in) :: nsoil_PEM ! # of soil layers in the PEM |
---|
[3149] | 38 | integer, intent(in) :: nsoil_PCM ! # of soil layers in the PCM |
---|
| 39 | real, dimension(ngrid,nsoil_PCM,nslope), intent(in) :: TI_PCM ! Thermal inertia in the PCM [SI] |
---|
[2794] | 40 | |
---|
[3076] | 41 | real, dimension(ngrid,nsoil_PEM,nslope), intent(inout) :: TI_PEM ! Thermal inertia in the PEM [SI] |
---|
| 42 | !======================================================================= |
---|
[2794] | 43 | ! local variables: |
---|
[3076] | 44 | integer :: ig, iloop, islope ! loop counters |
---|
[3206] | 45 | real :: alpha, lay1 ! coefficients for building layers |
---|
| 46 | real :: index_powerlaw ! coefficient to match the power low grid with the exponential one |
---|
[3076] | 47 | !======================================================================= |
---|
[2794] | 48 | ! 1. Depth coordinate |
---|
| 49 | ! ------------------- |
---|
[3123] | 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) |
---|
[3076] | 52 | lay1 = 2.e-4 |
---|
| 53 | alpha = 2 |
---|
[3149] | 54 | do iloop = 0,nsoil_PCM - 1 |
---|
[3123] | 55 | mlayer_PEM(iloop) = lay1*(1+iloop**2.9*(1-exp(-real(iloop)/20.))) |
---|
[3076] | 56 | enddo |
---|
[2794] | 57 | |
---|
[3123] | 58 | do iloop = 0, nsoil_PEM-1 |
---|
[3149] | 59 | if(lay1*(alpha**(iloop-0.5)) > mlayer_PEM(nsoil_PCM-1)) then |
---|
[3123] | 60 | index_powerlaw = iloop |
---|
| 61 | exit |
---|
| 62 | endif |
---|
[3076] | 63 | enddo |
---|
[2794] | 64 | |
---|
[3149] | 65 | do iloop = nsoil_PCM, nsoil_PEM-1 |
---|
| 66 | mlayer_PEM(iloop) = lay1*(alpha**(index_powerlaw + (iloop-nsoil_PCM)-0.5)) |
---|
[3123] | 67 | enddo |
---|
| 68 | |
---|
| 69 | ! Build layer_PEM() |
---|
| 70 | do iloop=1,nsoil_PEM-1 |
---|
| 71 | layer_PEM(iloop)=(mlayer_PEM(iloop)+mlayer_PEM(iloop-1))/2. |
---|
| 72 | enddo |
---|
| 73 | layer_PEM(nsoil_PEM)=2*mlayer_PEM(nsoil_PEM-1) - mlayer_PEM(nsoil_PEM-2) |
---|
| 74 | |
---|
| 75 | |
---|
[2794] | 76 | ! 2. Thermal inertia (note: it is declared in comsoil_h) |
---|
| 77 | ! ------------------ |
---|
[3076] | 78 | do ig = 1,ngrid |
---|
| 79 | do islope = 1,nslope |
---|
[3149] | 80 | do iloop = 1,nsoil_PCM |
---|
| 81 | TI_PEM(ig,iloop,islope) = TI_PCM(ig,iloop,islope) |
---|
[2794] | 82 | enddo |
---|
[3149] | 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) |
---|
[3076] | 86 | enddo |
---|
| 87 | endif |
---|
| 88 | enddo |
---|
| 89 | enddo |
---|
[2835] | 90 | |
---|
[2895] | 91 | ! 3. Index for subsurface layering |
---|
| 92 | ! ------------------ |
---|
[3076] | 93 | index_breccia = 1 |
---|
| 94 | do iloop = 1,nsoil_PEM - 1 |
---|
| 95 | if (depth_breccia >= layer_PEM(iloop)) then |
---|
| 96 | index_breccia = iloop |
---|
| 97 | else |
---|
| 98 | exit |
---|
| 99 | endif |
---|
| 100 | enddo |
---|
[2895] | 101 | |
---|
[3076] | 102 | index_bedrock = 1 |
---|
| 103 | do iloop = 1,nsoil_PEM - 1 |
---|
| 104 | if (depth_bedrock >= layer_PEM(iloop)) then |
---|
| 105 | index_bedrock = iloop |
---|
| 106 | else |
---|
| 107 | exit |
---|
| 108 | endif |
---|
| 109 | enddo |
---|
[2895] | 110 | |
---|
[3076] | 111 | END SUBROUTINE soil_settings_PEM |
---|
[2895] | 112 | |
---|
[3076] | 113 | END MODULE soil_settings_PEM_mod |
---|