Changeset 3991 for trunk/LMDZ.COMMON/libf/evolution/soil.F90
- Timestamp:
- Dec 16, 2025, 4:39:24 PM (4 weeks ago)
- File:
-
- 1 edited
-
trunk/LMDZ.COMMON/libf/evolution/soil.F90 (modified) (5 diffs)
Legend:
- Unmodified
- Added
- Removed
-
trunk/LMDZ.COMMON/libf/evolution/soil.F90
r3989 r3991 1 1 MODULE soil 2 3 implicit none 4 5 integer, parameter :: nsoilmx_PEM = 68 ! number of layers in the PEM 2 !----------------------------------------------------------------------- 3 ! NAME 4 ! soil 5 ! 6 ! DESCRIPTION 7 ! Soil state and properties for PEM. 8 ! 9 ! AUTHORS & DATE 10 ! L. Lange, 04/2023 11 ! JB Clement, 2023-2025 12 ! 13 ! NOTES 14 ! 15 !----------------------------------------------------------------------- 16 17 ! DECLARATION 18 ! ----------- 19 implicit none 20 21 ! MODULE PARAMETERS 22 ! ----------------- 23 integer, parameter :: nsoilmx_PEM = 68 ! number of layers in the PEM 24 25 ! MODULE VARIABLES 26 ! ----------------- 6 27 real, allocatable, dimension(:) :: layer_PEM ! soil layer depths [m] 7 28 real, allocatable, dimension(:) :: mlayer_PEM ! soil mid-layer depths [m] … … 25 46 logical :: reg_thprop_dependp ! thermal properites of the regolith vary with the surface pressure 26 47 27 !=======================================================================28 48 contains 29 !======================================================================= 30 49 !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ 50 51 !======================================================================= 31 52 SUBROUTINE ini_soil(ngrid,nslope) 32 33 implicit none 34 53 !----------------------------------------------------------------------- 54 ! NAME 55 ! ini_soil 56 ! 57 ! DESCRIPTION 58 ! Allocate soil arrays based on grid dimensions. 59 ! 60 ! AUTHORS & DATE 61 ! L. Lange, 2023 62 ! JB Clement, 2023-2025 63 ! 64 ! NOTES 65 ! 66 !----------------------------------------------------------------------- 67 68 ! DECLARATION 69 ! ----------- 70 implicit none 71 72 ! ARGUMENTS 73 ! --------- 35 74 integer, intent(in) :: ngrid ! number of atmospheric columns 36 75 integer, intent(in) :: nslope ! number of slope within a mesh 37 76 77 ! CODE 78 ! ---- 38 79 allocate(layer_PEM(nsoilmx_PEM)) 39 80 allocate(mlayer_PEM(0:nsoilmx_PEM - 1)) … … 49 90 50 91 END SUBROUTINE ini_soil 51 52 !======================================================================= 53 92 !======================================================================= 93 94 !======================================================================= 54 95 SUBROUTINE end_soil 55 56 implicit none 57 96 !----------------------------------------------------------------------- 97 ! NAME 98 ! end_soil 99 ! 100 ! DESCRIPTION 101 ! Deallocate all soil arrays. 102 ! 103 ! AUTHORS & DATE 104 ! L. Lange, 2023 105 ! JB Clement, 2023-2025 106 ! 107 ! NOTES 108 ! 109 !----------------------------------------------------------------------- 110 111 ! DECLARATION 112 ! ----------- 113 implicit none 114 115 ! CODE 116 ! ---- 58 117 if (allocated(layer_PEM)) deallocate(layer_PEM) 59 118 if (allocated(mlayer_PEM)) deallocate(mlayer_PEM) … … 69 128 70 129 END SUBROUTINE end_soil 130 !======================================================================= 71 131 72 132 !======================================================================= 73 133 SUBROUTINE set_soil(ngrid,nslope,nsoil_PEM,nsoil_PCM,TI_PCM,TI_PEM) 74 134 !----------------------------------------------------------------------- 135 ! NAME 136 ! set_soil 137 ! 138 ! DESCRIPTION 139 ! Initialize soil depths and properties. Builds layer depths using 140 ! exponential then power-law distribution; interpolates thermal inertia 141 ! from PCM to PEM grid. 142 ! 143 ! AUTHORS & DATE 144 ! E. Millour, 07/2006 145 ! L. Lange, 2023 146 ! JB Clement, 2023-2025 147 ! 148 ! NOTES 149 ! Modifications: 150 ! Aug.2010 EM: use NetCDF90 to load variables (enables using 151 ! r4 or r8 restarts independently of having compiled 152 ! the GCM in r4 or r8) 153 ! June 2013 TN: Possibility to read files with a time axis 154 ! The various actions and variable read/initialized are: 155 ! 1. Read/build layer (and midlayer) depth 156 ! 2. Interpolate thermal inertia and temperature on the new grid, if necessary 157 !----------------------------------------------------------------------- 158 159 ! DEPENDENCIES 160 ! ------------ 75 161 use iostart, only: inquire_field_ndims, get_var, get_field, inquire_field, inquire_dimension_length 76 162 77 implicit none 78 79 !======================================================================= 80 ! Author: LL, based on work by Ehouarn Millour (07/2006) 81 ! 82 ! Purpose: Read and/or initialise soil depths and properties 83 ! 84 ! Modifications: Aug.2010 EM: use NetCDF90 to load variables (enables using 85 ! r4 or r8 restarts independently of having compiled 86 ! the GCM in r4 or r8) 87 ! June 2013 TN: Possibility to read files with a time axis 88 ! 89 ! The various actions and variable read/initialized are: 90 ! 1. Read/build layer (and midlayer) depth 91 ! 2. Interpolate thermal inertia and temperature on the new grid, if necessary 92 !======================================================================= 93 94 !======================================================================= 95 ! arguments 96 ! --------- 97 ! inputs: 98 integer, intent(in) :: ngrid ! # of horizontal grid points 99 integer, intent(in) :: nslope ! # of subslope wihtin the mesh 100 integer, intent(in) :: nsoil_PEM ! # of soil layers in the PEM 101 integer, intent(in) :: nsoil_PCM ! # of soil layers in the PCM 102 real, dimension(ngrid,nsoil_PCM,nslope), intent(in) :: TI_PCM ! Thermal inertia in the PCM [SI] 103 104 real, dimension(ngrid,nsoil_PEM,nslope), intent(inout) :: TI_PEM ! Thermal inertia in the PEM [SI] 105 !======================================================================= 106 ! local variables: 163 ! DECLARATION 164 ! ----------- 165 implicit none 166 167 ! ARGUMENTS 168 ! --------- 169 integer, intent(in) :: ngrid ! # of horizontal grid points 170 integer, intent(in) :: nslope ! # of subslope wihtin the mesh 171 integer, intent(in) :: nsoil_PEM ! # of soil layers in the PEM 172 integer, intent(in) :: nsoil_PCM ! # of soil layers in the PCM 173 real, dimension(ngrid,nsoil_PCM,nslope), intent(in) :: TI_PCM ! Thermal inertia in the PCM [SI] 174 real, dimension(ngrid,nsoil_PEM,nslope), intent(inout) :: TI_PEM ! Thermal inertia in the PEM [SI] 175 176 ! LOCAL VARIABLES 177 ! --------------- 107 178 integer :: ig, iloop, islope ! loop counters 108 179 real :: alpha, lay1 ! coefficients for building layers 109 180 real :: index_powerlaw ! coefficient to match the power low grid with the exponential one 110 !======================================================================= 181 182 ! CODE 183 ! ---- 111 184 ! 1. Depth coordinate 112 185 ! ------------------- … … 173 246 174 247 END SUBROUTINE set_soil 248 !======================================================================= 175 249 176 250 END MODULE soil
Note: See TracChangeset
for help on using the changeset viewer.
