source: trunk/LMDZ.MARS/libf/phymars/soilwater.F90 @ 3225

Last change on this file since 3225 was 3223, checked in by llange, 12 months ago

Mars PCM
Adsorption: commit the last version developed by Pierre Yves Meslin
LL

File size: 96.5 KB
Line 
1subroutine soilwater(ngrid, nlayer, nq, nsoil, nqsoil, ptsrf, ptsoil, ptimestep,  &
2      exchange, qsat_surf, pq, pa, pb, pc, pd, pdqsdifpot, pqsurf,  &
3      pqsoil, pplev, rhoatmo, writeoutput, zdqsdifrego, zq1temp2)
4
5
6      use comsoil_h, only: igcm_h2o_vap_soil, igcm_h2o_ice_soil, igcm_h2o_vap_ads, layer, mlayer, choice_ads
7      use comcstfi_h
8      use tracer_mod
9      use surfdat_h, only: watercaptag ! use mis par AP15 essai
10      use geometry_mod, only: cell_area, latitude_deg
11      use write_output_mod, only: write_output
12implicit none
13
14! =================================================================================================  =
15! Description:
16!
17! This subroutine calculates the profile of water vapor, adsorbed water and ice in the regolith,
18! down to depth(nsoil) ~ 18m. It calculates the flux of water vapor (zdqsdifrego) between the
19! suburface and the atmosphere.
20!
21! Water vapor and nsoiladsorbed water are treated as two separate subsurface tracers that are connected
22! by adsorption / desorption coefficients, including an adsorption / desorption kinetic factor.
23!
24! The water adsorption isotherm is linear, with a saturation plateau corresponding to one monolayer
25! of adsorbed water molecules (i.e. an approximation of a Langmuir - type isotherm). The slope of the
26! isotherm is defined by the enthalpy of adsorption (in kJ / mol).
27! New since 2020 : the adsorption isotherm is now linear by part, to better simulate Type II or Type III isotherms (and not only Type I)
28!
29! The linearized adsorption - diffusion equation is solved first, then the water vapor pressure is
30! compared to the saturation vapor pressure, and if the latter is reached, ice is formed, the water vapor
31! pressure is fixed to its saturation value and the adsorption - diffusion equation is solved again.
32! If ice was already present, its sublimation or growth is calculated.
33!
34! This subroutine is called by vdifc.F if the parameters "regolithe" and "water" are set to true in
35! callphys.def. It returns the flux of water vapor that is injected into or out of the regolith,
36! and which serves as a boundary condition to calculate the vertical atmospheric profile of water vapor
37! in vdifc. It also calculates the exchange between the subsurface and surface ice.
38!
39! It requires three subsurface tracers to be defined in traceur.def :
40!    -  h2o_vap_soil (subsurface water vapor)
41!    -  h2o_ice_soil (subsurface ice)
42!    -  h2o_ads_vap  (adsorbed water)
43!
44! Author : Pierre - Yves Meslin (pmeslin@irap.omp.eu)
45!
46! =================================================================================================  =
47
48! Libraries :
49! ===========  =
50!#include "dimensions.h"
51!#include "dimphys.h"
52!#include "comsoil.h"
53#include "callkeys.h"
54!#include "comcstfi.h"
55!#include "tracer.h"
56!#include "watercap.h"
57
58
59! Arguments :
60! ======================================================================
61
62! Inputs :
63! ----------------------------------------------------------------------
64integer, intent(in) :: ngrid                    ! number of points in the model (all lat and long point in one 1D array)   
65integer, intent(in) :: nlayer                   ! number of layers
66integer, intent(in) :: nq                       ! number of tracers
67integer, intent(in) :: nsoil
68integer, intent(in) :: nqsoil
69logical, save :: firstcall_soil = .true.             ! triggers the initialization
70real, intent(in) :: ptsoil(ngrid, nsoil)         ! Subsurface temperatures
71real, intent(in) :: ptsrf(ngrid)                ! Surface temperature
72real, intent(in) :: ptimestep                   ! length of the timestep (unit depends on run.def file)
73logical, intent(in) :: exchange(ngrid)          ! logical :: array describing whether there is exchange with the atmosphere at the current timestep
74
75real, intent(in) :: qsat_surf(ngrid)           ! Saturation vapor pressure at the current surface temperature (formerly qsat)
76real, intent(in) :: pq(ngrid, nlayer, nq)       ! Tracer atmospheric mixing ratio
77real, intent(in) :: pa(ngrid, nlayer)           ! Coefficients za
78real, intent(in) :: pb(ngrid, nlayer)           ! Coefficients zb
79real, intent(in) :: pc(ngrid, nlayer)           ! Coefficients zc
80real, intent(in) :: pd(ngrid, nlayer)           ! Coefficients zd
81real, intent(in) :: pdqsdifpot(ngrid, nq)       ! Potential pdqsdif (without subsurface - atmosphere exchange)
82
83real, intent(in) :: pplev(ngrid, nlayer+1)      ! Atmospheric pressure levels
84real, intent(in) :: rhoatmo(ngrid)              ! Atmospheric air density (1st layer) (not used right now)
85logical, intent(in) :: writeoutput              ! just to check we are at the last subtimestep and we
86
87! Variables modified :
88! ----------------------------------------------------------------------
89real, intent(inout) :: pqsoil(ngrid, nsoil, 3*nq) ! Subsurface tracers (water vapor and ice)
90real, intent(in) :: pqsurf(ngrid)                 ! Water ice on the surface
91                                                  ! (in kg.m - 3 : m - 3 of pore air for water vapor and m - 3 of regolith for water ice and adsorbed water)
92! Outputs :
93! ----------------------------------------------------------------------
94real, intent(out) :: zdqsdifrego(ngrid)         ! Flux from subsurface (positive pointing outwards)
95real, intent(out) :: zq1temp2(ngrid)            ! Temporary atmospheric mixing ratio after exchange with subsurface (kg / kg)
96
97! Outputs for the output files
98real*8 :: preduite(ngrid, nsoil)                ! how close the water vapor density is to adsorption saturation
99integer :: exch(ngrid)                          ! translates the logical variable exchange into a -1 or 1
100real*8 :: mass_h2o(ngrid)                       ! Mass of subsurface water column at each point (kg.m - 2) (formerly mh2otot)
101real*8 :: mass_ice(ngrid)                       ! Mass of subsurface ice at each point (kg.m - 2) (formerly micetot)
102real*8 :: mass_h2o_tot                          ! Mass of subsurface water over the whole planet
103real*8 :: mass_ice_tot                          ! Mass of subsurface ice over the whole planet
104real*8 :: nsurf(ngrid)                          ! surface tracer density in kg/m^3
105real*8 :: close_out(ngrid, nsoil)           ! output for close_top and close_bottom
106
107! Local (saved) variables :
108! ======================================================================
109
110real*8 :: P_sat_soil(ngrid, nsoil)              ! water saturation pressure of subsurface cells (at miD-layers) (formerly Psatsoil)
111real*8 :: nsatsoil(ngrid, nsoil)                ! gas number density at saturation pressure
112real*8, allocatable, save :: znsoil(:, :)       ! Water vapor soil concentration (kg.m - 3 of pore air)
113real*8 :: znsoilprev(ngrid, nsoil)              ! value of znsoil in the previous timestep
114real*8 :: znsoilprev2(ngrid, nsoil)             ! second variable for the value of znsoil in the previous timestep
115real*8 :: zdqsoil(ngrid, nsoil)                 ! Increase of pqsoil if sublimation of ice
116real*8, allocatable, save :: ice(:, :)          ! Ice content of subsurface cells (at miD-layers) (kg / m3)
117real*8 :: iceprev(ngrid, nsoil)
118logical :: ice_index_flag(nsoil)                ! flag for ice presence
119real*8, allocatable, save :: adswater(:, :)     ! Adsorbed water in subsurface cells (at miD-layers) (...)
120real*8 :: adswater_temp(ngrid, nsoil)           ! temprory variable for adsorbed water
121logical, allocatable, save  :: over_mono_sat_flag(:, :) ! Formerly ads_water_saturation_flag_2(nsoil) (see descritption of the variable recompute_cell_ads_flag for an explanation ! ARNAU
122real*8 :: adswprev(ngrid, nsoil)
123logical :: recompute_cell_ads_flag(nsoil) ! ARNAU
124! Formerly ads_water_saturation_flag_1 but with a twist. This variable now
125! checks whether coefficients have changed and need recomputing. If the cell
126! is seen to be over the monolayer saturation level (i.e. the cell fulfills the
127! condition adswater_temp(ig, ik) > adswater_sat_mono) but the coefficients
128! have been computed assuming that the cell is below the monolayer saturation
129! layer (i.e. the variable over_mono_sat_flag(ig, ik) had been set to .false.),
130! then the cell needs to have its coefficients recomputed according to the
131! previous condition: adswater_temp(ig, ik) > adswater_sat_mono. Then,
132! the variable recompute_cell_ads_flag(ik) becomes .true.. ! ARNAU
133real*8, save :: adswater_sat_mono                    ! Adsorption monolayer saturation value (kg.m - 3 of regolith) ! ARNAU
134real*8 :: delta0(ngrid)                         ! Coefficient delta(0) modified
135real*8 :: alpha0(ngrid)
136real*8 :: beta0(ngrid)
137
138! Adsorbtion/Desorption variables and parameters
139real*8 :: Ka(ngrid, nsoil)            ! Adsorption time constant (s - 1) before monolayer saturation (first segment of the bilinear function)
140real*8 :: Kd(ngrid, nsoil)            ! Desorption time constant (s - 1) before monolayer saturation (first segment of the bilinear function)
141real*8 :: k_ads_eq(ngrid, nsoil)      ! Equilibrium adsorption coefficient (formerly kads) before monolayer saturation (first segment of the bilinear function); unitless (converts kg/m3 into kg/m3)
142real*8 :: Ka2(ngrid, nsoil)           ! Adsorption time constant (s - 1) after monolayer saturation (second segment of the bilinear function) ! modified 2020
143real*8 :: Kd2(ngrid, nsoil)           ! Desorption time constant (s - 1) after monolayer saturation (second segment of the bilinear function) ! modified 2020
144real*8 :: k_ads_eq2(ngrid, nsoil)     ! Equilibrium adsorption coefficient (formerly kads) after monolayer saturation (second segment of the bilinear function); unitless ! modified 2020
145real*8 :: c0(ngrid, nsoil)            ! Intercept of the second line in the bilinear function ! modified 2020
146real*8, parameter :: kinetic_factor = 0.01      ! (inverse of) Characteristic time to reach adsorption equilibrium (s - 1):
147                                                ! fixed arbitrarily when kinetics factors are unknown
148                                                ! possible values: ! = 1 / 1800 s - 1 ! / 1.16D-6 /  !( =  10 earth days)! / 1.D0 /  ! / 1.2D-5 /  !
149                                         
150real*8, allocatable, save :: ztot1(:, :)  ! Total (water vapor +  ice) content (kg.m - 3 of soil)
151real*8 :: dztot1(nsoil)
152real*8 :: h2otot(ngrid, nsoil)      ! Total water content (ice +  water vapor +  adsorbed water)
153real*8 :: flux(ngrid, 0:nsoil - 1)  ! Fluxes at interfaces (kg.m - 2.s - 1) (positive = upward)
154real*8 :: rho(ngrid)                ! Air density (first subsurface layer)
155real*8 :: rhosurf(ngrid)            ! Surface air density
156
157
158! Porosity and tortuosity variables
159real*8, allocatable, save :: porosity_ice_free(:, :)  ! Porosity with no ice present (formerly eps0)
160real*8, allocatable, save :: porosity(:, :)           ! Porosity (formerly eps)
161real*8 :: porosity_prev(ngrid, nsoil)                 ! Porosity from previous timestep
162real*8, allocatable, save :: tortuosity(:, :)         ! Tortuosity factor (formerly tortuo)
163
164real*8 :: saturation_water_ice(ngrid, nsoil)          ! Water pore ice saturation level (formerly Sw)
165real*8 :: saturation_water_ice_inter(ngrid, nsoil)          ! Water pore ice saturation level at the interlayer
166
167! Diffussion coefficients
168real*8 :: D_mid(ngrid, nsoil)       ! Midlayer diffusion coefficients
169real*8 :: D_inter(ngrid, nsoil)     ! Interlayer diffusion coefficients (formerly D)
170real*8, allocatable, save :: D0(:, :)     ! Diffusion coefficient prefactor :
171                                          ! If (medium = 1) : D0 = value of D_mid for saturation_water_ice = 0, ie. poro / tortuo * Dm (in m2 / s)
172                                          ! If (medium = 2) : D0 = porosity_tortuosity factor (dimensionless)
173real*8 :: omega(ngrid, nsoil)       ! H2O - CO2 collision integral
174real*8 :: vth(ngrid, nsoil)         ! H2O thermal speed
175real*8, parameter :: Dk0 = 0.459D0  ! Knudsen diffusion coefficient (for saturation_water_ice = 0) Value for a (5 / 1) bidispersed random assembly of spheres
176                                    ! (dimensionless, ie. divided by thermal speed and characteristic meshsize of the porous medium)
177real*8 :: Dk(ngrid, nsoil)          ! Knudsen diffusion coefficient (divided by porosity / tortuosity factor)
178real*8 :: Dk_inter(ngrid, nsoil)          ! Knudsen diffusion coefficient at the interlayer
179real*8 :: Dm(ngrid, nsoil)          ! Molecular diffusion coefficient
180real*8 :: Dm_inter(ngrid, nsoil)          ! Molecular diffusion coefficient at the interlayer
181
182real*8, parameter :: choke_fraction = 0.8D0  ! fraction of ice that prevents further diffusion
183logical, allocatable, save :: close_top(:, :)         ! closing diffusion at the top of a layer if ice rises over saturation
184logical, allocatable, save :: close_bottom(:, :)      ! closing diffusion at the bottom of a layer if ice risies over saturation
185logical, parameter :: print_closure_warnings = .true.
186
187! Coefficients for the Diffusion calculations
188real*8 :: A(ngrid, nsoil)           ! Coefficient in the diffusion formula
189real*8 :: B(ngrid, nsoil)           ! Coefficient in the diffusion formula
190real*8 :: C(ngrid, nsoil)           ! Coefficient in the diffusion formula
191real*8 :: E(ngrid, nsoil)           ! Coefficient in the diffusion formula (before monolayer saturation) ! added 2020
192real*8 :: F(ngrid, nsoil)           ! Coefficient in the diffusion formula (before monolayer saturation) ! added 2020
193real*8 :: E2(ngrid, nsoil)           ! Coefficient in the diffusion formula (after monolayer saturation) ! added 2020
194real*8 :: F2(ngrid, nsoil)           ! Coefficient in the diffusion formula (after monolayer saturation) ! added 2020
195real*8, allocatable, save :: zalpha(:, :) ! Alpha coefficients
196real*8, allocatable, save :: zbeta(:, :)  ! Beta coefficients
197real*8 :: zdelta(ngrid, nsoil - 1)        ! Delta coefficients
198
199! Distances between layers
200real*8, allocatable, save :: interlayer_dz(:, :)      ! Distance between the interlayer points in m (formerly interdz)
201real*8, allocatable, save :: midlayer_dz(:, :)        ! Distance between the midcell points in m (formerly middz)
202
203real*8 :: nsat(ngrid, nsoil)                    ! amount of water vapor at (adsorption) monolayer saturation ! modified 2020
204
205real*8, allocatable, save :: meshsize(:, :)     ! scaling/dimension factor for the por size
206real*8, allocatable, save :: rho_soil(:, :)     ! Soil density (bulk -  kg / m3) (formerly rhosoil)
207real*8, allocatable, save :: cste_ads(:, :)     ! Prefactor of adsorption coeff. k
208
209! general constants
210real*8, parameter :: kb = 1.38065D-23     ! Boltzmann constant
211real*8, parameter :: mw = 2.988D-26             ! Water molecular weight
212!real*8, parameter :: rho_H2O_ice = 920.D0    ! Ice density
213
214! adsorption coefficients
215real*8, parameter :: enthalpy_ads = 35.D3 ! Enthalpy of adsorption (J.mole - 1) options 21.D3, 35.D3, and 60.D3
216real*8, parameter :: enthalpy_ads2 = 21.D3 ! Enthalpy of adsorption (J.mole - 1) options 21.D3, 35.D3, and 60.D3 for the second segment ! added 2020
217real*8, parameter :: DeltaQ_ads = 21.D3 ! 10.D3 ! This is the DeltaQ_ads = heat of adsorption - enthalpy of liquefaction/sublimation = E1 - EL and may be equal to RT*ln(c), where c is the BET constant (from BET1938). This is for the first segment (J.mole - 1) ! added 2020
218real*8, parameter :: DeltaQ_ads2 = 21.D3 ! 10.D3 ! This is the DeltaQ_ads = heat of adsorption - enthalpy of liquefaction/sublimation = E1 - EL and may be equal to RT*ln(c), where c is the BET constant (from BET1938). This is for the second segment (J.mole - 1) ! added 2020
219real*8, parameter :: tau0 = 1D-14
220real*8, parameter :: S = 17.D3            ! Soil specific surface area (m2.kg - 1 of solid) options: 17.D3 and 106.D3
221real*8, parameter:: Sm = 10.6D-20         ! Surface of the water molecule (m2) (only needed in the theoretical formula which is not used right now)
222
223
224! Reference values for choice_ads = 2
225real*8, parameter :: Tref = 243.15D0
226real*8, parameter :: nref = 2.31505631D-6 ! calculated as 18.D-3 / (8.314D0 * Tref) * 0.26D0 ! not used anymore (for the time being)
227real*8, parameter :: Kd_ref = 0.0206D0   ! Not used for the time being (would require specific measurements to be known, but it is rarely measured)
228real*8, parameter :: Ka_ref = 2.4D-4     ! Not used for the time being
229! real*8, parameter :: Kref = 6.27D6 ! Value from data analysis from Pommerol2009 ! Old value 1.23D7        ! Volcanic tuff @ 243.15 K (obtained at low P / Psat) ! ARNAU
230! real*8, parameter :: Kref2 = 5.95D-7 ! Value from data analysis Pommerol2009 ! ARNAU
231real*8, parameter :: Kref = 0.205D-6 ! Value obtained from the fit of all adsorption data from Pommmerol (2009) (see Arnau's report - p.28: = yp/xp = 2.6998E-7/3.5562E-2, divided by psat(243.15K=37 Pa; will need to be multiplied by RT/M to be unitless before multiplying by znsoil, which in kg(water)/m3(air) and not in pascals)
232real*8, parameter :: Kref2 = 0.108D-7 ! Value obtained from the fit of all adsorption data from Pommmerol (2009) (see Arnau's report - p.28 = m2; divided by psat(243.15K)=37 Pa)
233
234
235logical :: recompute_all_cells_ads_flag ! Old saturation_flag but with a behaviour change ! ARNAU
236! The old saturation_flag was used to check whether the cell was saturated or
237! not, and therefore to assign the saturation value adswater(ig, ik) =
238! adswater_sat instead of the usual adswater(ig, ik) = adswater_temp(ig, ik)
239! The old routine using saturation_flag did not require that the A, B, C,...
240! adsorption coefficients be computed, but the new soilwater subroutine
241! does. Therefore, the variable recompute_all_cells_ads_flag checks whether
242! there is a cell of a column that requires recomputing. If at least one cell
243! requires recomputing (i.e. recompute_cell_ads_flag(ik) is .true.), then
244! recompute_all_cells_ads_flag becomes true, and the adsorption coefficients,
245! as well as the recursive equation (appearing in this code as `backward
246! loop over all layers`) ! ARNAU
247logical :: sublimation_flag
248logical :: condensation_flag(nsoil)
249
250! variables and parameters for the H2O map
251integer, parameter :: n_long_H2O = 66!180          ! number of longitudes of the new H2O map
252integer, parameter :: n_lat_H2O = 50 !87            ! number of latitudes of the new H2O map
253real*8, parameter :: rho_H2O_ice = 920.0D0      ! Ice density (formerly rhoice)
254real :: latH2O(n_lat_H2O*n_long_H2O)            ! Latitude at H2O map points
255real :: longH2O(n_lat_H2O*n_long_H2O)           ! Longitude at H2O map points
256real :: H2O_depth_data(n_lat_H2O*n_long_H2O)    ! depth of the ice layer in in g/cm^2 at H2O map points
257real :: H2O_cover_data(n_lat_H2O*n_long_H2O)    ! H2O content of the cover layer at H2O map points (not used yet)
258real :: dataH2O(n_lat_H2O*n_long_H2O)           ! H2O content of the ice layer at H2O map points
259real :: latH2O_temp(n_lat_H2O*n_long_H2O)       ! intermediate variable
260real :: longH2O_temp(n_lat_H2O*n_long_H2O)      ! intermediate variable
261real :: dataH2O_temp(n_lat_H2O*n_long_H2O)      ! intermediate variable
262real :: H2O_depth_data_temp(n_lat_H2O*n_long_H2O)     ! intermediate variable
263real, allocatable, save :: H2O(:)                            ! H2O map interpolated at GCM grid points (in wt%)
264real, allocatable, save :: H2O_depth(:)                      ! H2O map ice depth interpolated at GCM in g/cm^2
265
266! variables for the 1D case
267real*8, parameter :: mmr_h2o = 0.6D-4     ! Water vapor mass mixing ratio (for initialization only)
268
269real*8 :: diff(ngrid, nsoil)              ! difference between total water content and ice + vapor content
270                                          ! (only used for output, inconstistent: should just be adswater)
271real :: var_flux_soil(ngrid, nsoil)
272
273! Loop variables and counters
274integer :: ig, ik, i, j                   ! loop variables
275logical :: output_trigger                 ! used to limit amount of written output
276integer, save :: n                        ! number of runs of this subroutine
277integer :: sublim_n                       ! number of sublimation loop runs
278integer :: satur_mono_n                   ! number of monolayer saturation loop runs ! added 2020
279
280
281if (.not. allocated(znsoil)) then
282      allocate( znsoil(ngrid, nsoil) )
283      allocate( ice(ngrid, nsoil) )
284      allocate( adswater(ngrid, nsoil) )
285      allocate( ztot1(ngrid, nsoil) )
286      allocate( porosity_ice_free(ngrid, nsoil) )
287      allocate( porosity(ngrid, nsoil) )
288      allocate( tortuosity(ngrid, nsoil) )
289      allocate( D0(ngrid, nsoil) )
290      allocate( interlayer_dz(ngrid, nsoil) )
291      allocate( midlayer_dz(ngrid, 0:nsoil) )
292!      allocate( zalpha(ngrid, nsoil-1) )
293!      allocate( zbeta(ngrid, nsoil-1) )
294      allocate( zalpha(ngrid, nsoil) )    ! extra entry to match output format
295      allocate( zbeta(ngrid, nsoil) )     ! extra entry to match output format
296      allocate( meshsize(ngrid, nsoil) )
297      allocate( rho_soil(ngrid, nsoil) )
298      allocate( cste_ads(ngrid, nsoil) )
299      allocate( H2O(ngrid) )
300      allocate( H2O_depth(ngrid) )
301      allocate( close_top(ngrid, nsoil) )
302      allocate( close_bottom(ngrid, nsoil) )
303      allocate( over_mono_sat_flag(ngrid, nsoil) )
304endif
305
306! Used commons
307! mugaz ! Molar mass of the atmosphere (g.mol-1) ~43.49
308
309
310
311! Initialisation :
312! =================
313
314if (firstcall_soil) then
315      n = 0
316      close_top = .false.
317      close_bottom = .false.
318      print *, "adsorption enthalpy first segment: ", enthalpy_ads
319      print *, "adsorption enthalpy second segment: ", enthalpy_ads2
320      print *, "adsorption BET constant C first segment: ", DeltaQ_ads
321      print *, "adsorption BET constant C second segment: ", DeltaQ_ads2
322      print *, "specific surface area: ", S
323
324      do ig = 1, ngrid
325!            if(.not.watercaptag(ig)) then
326
327            ! Initialization of soil parameters
328            ! ===================================
329
330                  ! initialize the midlayer distances
331                  midlayer_dz(ig, 0) = mlayer(0)
332                 
333                  do ik = 1, nsoil - 1
334                        midlayer_dz(ig, ik) = mlayer(ik) - mlayer(ik - 1)
335                  enddo
336
337                  ! initialize the interlayer distances
338                  interlayer_dz(ig, 1) = layer(1)
339                  do ik = 2, nsoil
340                        interlayer_dz(ig, ik) = layer(ik) - layer(ik - 1)
341                  enddo
342
343                  ! Choice of porous medium and D0
344                  ! ===============================  =
345                  do ik = 1, nsoil
346                       
347                        ! These properties are defined here in order to enable custom profiles
348                        porosity_ice_free(ig, ik) = 0.45D0
349                        tortuosity(ig, ik) = 1.5D0
350                        rho_soil(ig, ik) = 1.3D3 ! in kg/m3 of regolith (incl. porosity)
351
352                        meshsize(ig, ik) = 5.D-6  ! Characteristic size 1e - 6 = 1 micron : diameter(mode 5) = 10 microns, diameter(mode 1) = 1 microns
353                                                  ! Example : with meshsize = 5.E - 6 : grain sizes = 50 and 10 microns
354                        D0(ig, ik) = porosity_ice_free(ig, ik) / tortuosity(ig, ik)
355
356                  enddo
357
358                  ! Choice of adsorption coefficients
359                  ! ===================================
360                    adswater_sat_mono = 2.6998D-7 * S * rho_soil(ig,1)  ! Unit = kg/m3; From best fit of all adsoprtion data from Pommerol et al. (2009) - See Arnau's report                   
361                  ! adswater_sat_mono = 0.8D-2 * S / 13.7D3 * rho_soil(ig, 1)  ! Unit = kg/m3 ; Experimental value for volcanic tuff (Pommerol et al., 2009)
362                  ! theoretical formula is = rho_soil(ig, 1) * S / Sm * mw
363
364                  ! Numerical values above are for SSA = 13.7 m2 / g and plateau = 0.8wt%
365                  !                       adswater_sat_mono = 6D-2 * rho_soil(ig, 1)            ! Experimental value for JSC Mars - 1 (Pommerol et al., 2009)
366                  ! with SSA = 106 m2 / g and plateau (net) = 6wt%
367
368                  ! Initialisation of ice, water vapor and adsorbed water profiles
369                  ! ===============================================================  =
370                  do ik = 1, nsoil
371                        ! Initialisation of pqsoil(t = 0)
372                        ! in 1D simulations the initialization happens here
373!                       if (ngrid.eq.1) then
374!                             znsoil(ig, ik) = mmr_h2o * mugaz * 1.D-3 * pplev(ig, 1) &
375!                              / (8.314D0 * ptsoil(ig, nsoil - 4))
376!                              !                                   znsoil(ig, ik) = 0.
377!                              ice(ig, ik) = 0.D0  ! 0.1D0 !414.D0
378                             
379!                              saturation_water_ice(ig, ik) = min(ice(ig, ik) / (rho_H2O_ice * porosity_ice_free(ig, ik)), 0.999D-0)
380!                              porosity(ig, ik) = porosity_ice_free(ig, ik) * (1.D0 - saturation_water_ice(ig, ik))
381                             
382!                              if (choice_ads.eq.1) then
383!                                    vth(ig, ik) = dsqrt(8.D0 * 8.314D0 * ptsoil(ig, nsoil - 4) &
384!                                          / (pi * 18.D-3))
385
386                                    ! first segment of the bilinear function
387!                                    k_ads_eq(ig, ik) = rho_soil(ig, ik) * S * vth(ig, ik) / 4.D0 &
388!                                          / (dexp(-enthalpy_ads &
389!                                          / (8.314D0 * ptsoil(ig, nsoil - 4))) / tau0)
390!                                    Kd(ig, ik) = kinetic_factor  &
391!                                          / (1.D0 + k_ads_eq(ig, ik) / porosity(ig, ik))
392!                                    Ka(ig, ik) = kinetic_factor * k_ads_eq(ig, ik) /  &
393!                                          (1.D0 + k_ads_eq(ig, ik) / porosity(ig, ik))
394
395                                    ! second segment of the bilinear function ! added 2020
396!                                    k_ads_eq2(ig, ik) = rho_soil(ig, ik) * S * vth(ig, ik) / 4.D0 &
397!                                          / (dexp(-enthalpy_ads2 &
398!                                          / (8.314D0 * ptsoil(ig, nsoil - 4))) / tau0)
399!                                    Kd2(ig, ik) = kinetic_factor  &
400!                                          / (1.D0 + k_ads_eq2(ig, ik) / porosity(ig, ik))
401!                                    Ka2(ig, ik) = kinetic_factor * k_ads_eq2(ig, ik) /  &
402!                                          (1.D0 + k_ads_eq2(ig, ik) / porosity(ig, ik))
403!
404!                              elseif (choice_ads.eq.2) then  ! par rapport \E0 valeur experimentale de reference
405!                                    !                                         Kd(ig, ik) = Kd_ref * exp(-enthalpy_ads / 8.314 *
406!                                    !     &                                         (1. / ptsoil(ig, nsoil - 4) - 1. / Tref))
407!                                    !                                         Ka(ig, ik) = rho_soil(ig, ik) * Ka_ref *
408!                                    !     &                                         sqrt(ptsoil(ig, nsoil - 4) / Tref) / nref!
409!
410!                                    ! first segment of the bilinear function 
411!                                    k_ads_eq(ig, ik) = 8.314D0 * ptsoil(ig,ik)/18.D-3 * Kref * dexp(DeltaQ_ads / 8.314D0 * (1.D0 / ptsoil(ig, ik) - 1.D0 / Tref)) ! Changed enthalpy_ads to DeltaQ_ads to ensure correct dependance/behaviour of k_ads_eq with temperature ; prefactor RT/M is to convert Kref in proper units to use with znsoil
412!                                         
413!                                    Kd(ig, ik) = kinetic_factor / (1.D0 + k_ads_eq(ig, ik) / porosity(ig, ik))
414!                                         
415!                                    Ka(ig, ik) = kinetic_factor * k_ads_eq(ig, ik) / (1.D0 + k_ads_eq(ig, ik) / porosity(ig, ik))
416!
417!                                    ! second segment of the bilinear function ! added 2020
418!                                    k_ads_eq2(ig, ik) = 8.314D0 * ptsoil(ig,ik)/18.D-3 * Kref2 * dexp(DeltaQ_ads2 / 8.314D0 * (1.D0 / ptsoil(ig, ik) - 1.D0 / Tref)) ! Changed enthalpy_ads2 to DeltaQ_ads2 to ensure correct dependance/behaviour of k_ads_eq with temperature
419!                                         
420!                                    Kd2(ig, ik) = kinetic_factor / (1.D0 + k_ads_eq2(ig, ik) / porosity(ig, ik))
421!                                         
422!                                    Ka2(ig, ik) = kinetic_factor * k_ads_eq2(ig, ik) / (1.D0 + k_ads_eq2(ig, ik) / porosity(ig, ik))
423!                              endif
424                               
425!                              if (k_ads_eq(ig,ik) * znsoil(ig, ik) .gt. adswater_sat_mono) then ! modified 2024, after correction of c0 expression
426!                                 c0(ig, ik) = ( 1.D0 - (k_ads_eq2(ig, ik) / k_ads_eq(ig, ik))) * adswater_sat_mono
427!                                 adswater(ig, ik) = c0(ig,ik) + k_ads_eq2(ig, ik) * znsoil(ig, ik)
428!                              else
429!                                 adswater(ig, ik) = k_ads_eq(ig, ik) * znsoil(ig, ik)
430!                              endif
431                             
432!                        else  ! in 3D simulations initialisation happens with newstart.F
433                              znsoil(ig, ik) = pqsoil(ig, ik, igcm_h2o_vap_soil)
434                              ice(ig, ik) = pqsoil(ig, ik, igcm_h2o_ice_soil)
435                              adswater(ig, ik) = pqsoil(ig, ik, igcm_h2o_vap_ads)
436!                        endif
437
438                        saturation_water_ice(ig, ik) = min(ice(ig, ik) / (rho_H2O_ice * porosity_ice_free(ig, ik)), 0.999D-0)
439                       
440                        if (saturation_water_ice(ig, ik).lt.0.D0) then
441                              print *, "WARNING!! soil water ice negative at ", ig, ik
442                              print *, "saturation value: ", saturation_water_ice(ig, ik)
443                              print *, "setting saturation to 0"
444                              saturation_water_ice(ig, ik) = max(saturation_water_ice(ig, ik), 0.D0)
445                        endif
446                       
447                        porosity(ig, ik) = porosity_ice_free(ig, ik) * (1.D0 - saturation_water_ice(ig, ik))
448                                               
449                        ! Initialise soil as if all points where below the monolayer saturation level (added 2020) => now depends on value of adswater (modified in 2024)
450                        if (adswater(ig,ik).gt.adswater_sat_mono) then
451                           over_mono_sat_flag(ig, ik) = .true. ! 
452                        else
453                            over_mono_sat_flag(ig, ik) = .false.   
454                        endif     
455
456                  enddo
457!            endif
458      enddo
459     
460      print *, "initializing H2O data"
461!     
462!      ! 1.6 intitalization of the water content
463!      open(40,file='H2O_data')
464!      do i = 1, n_lat_H2O*n_long_H2O
465!            read(40,*) latH2O_temp(i), longH2O_temp(i), H2O_cover_data(i), dataH2O_temp(i), H2O_depth_data_temp(i)
466!            ! write(*, *) 'depth data ', H2O_depth_data(i)
467!      enddo
468!      close(40) 
469!     
470!      ! 1.6.2. put latH2O, longH2O and dataH2O in the right format to pass on to mapTh
471!      ! in the datafile the latitudes are listed from negative to positive, but mapTh needs them the other way around
472!      do i = 0, n_lat_H2O - 1
473!            do j = 1, n_long_H2O
474!                  latH2O( n_long_H2O * i + j ) = latH2O_temp( n_long_H2O * (n_lat_H2O - (i + 1)) + j)
475!                  longH2O( n_long_H2O * i + j ) = longH2O_temp( n_long_H2O * (n_lat_H2O - (i + 1)) + j)
476!                  dataH2O( n_long_H2O * i + j ) = dataH2O_temp( n_long_H2O * (n_lat_H2O - (i + 1)) + j)
477!                  H2O_depth_data( n_long_H2O * i + j ) = H2O_depth_data_temp( n_long_H2O * (n_lat_H2O - (i + 1)) + j)
478!            enddo
479!      enddo
480!     
481!      call mapTh(latH2O, longH2O, dataH2O, n_long_H2O, n_lat_H2O, ngrid, H2O)
482
483!      call mapTh(latH2O, longH2O, H2O_depth_data, n_long_H2O, n_lat_H2O, ngrid, H2O_depth)
484!     
485      do ig = 1, ngrid
486            ! convert depth from g/cm^2 to m
487!            print *, H2O_depth(ig), rho_soil(ig, 1)
488!            H2O_depth(ig) = H2O_depth(ig) * 10 / rho_soil(ig, 1)
489!            if ( (latitude_deg(ig).lt.40).and.(latitude_deg(ig).gt.-40)) then
490!                  H2O(ig) = 0.
491!            endif
492           
493            output_trigger = .true.
494            do ik = 1, nsoil
495!                  if (H2O_depth(ig).le.layer(ik)) then
496!                        ice(ig, ik) = min(H2O(ig), rho_H2O_ice * porosity_ice_free(ig, ik))
497!                        if (output_trigger) then
498!                              print *, "ice set at: ", ig, ik, "to :", ice(ig,ik), "depth: ", H2O_depth(ig)
499!                              output_trigger = .false.
500!                        endif
501!                  endif
502                 
503                  saturation_water_ice(ig, ik) = min(ice(ig, ik) / (rho_H2O_ice * porosity_ice_free(ig, ik)), 0.999D-0)
504                  saturation_water_ice(ig, ik) = max(saturation_water_ice(ig, ik), 0.D0)
505                  porosity(ig, ik) = porosity_ice_free(ig, ik) * (1.D0 - saturation_water_ice(ig, ik))
506            enddo
507      enddo
508     
509      print *, "Kinetic factor: ", kinetic_factor
510      if (kinetic_factor.eq.0) then
511            print *, "WARNING: adsorption is switched off"
512      endif
513      firstcall_soil = .false.
514endif  ! of "if firstcall_soil"
515
516
517! Computation of new (new step) transport coefficients :
518! =======================================================  =
519
520do ig = 1, ngrid  ! loop over all gridpoints
521      if(.not.watercaptag(ig)) then  ! if there is no polar cap
522
523            do ik = 1, nsoil
524
525                  ! calculate Water saturation level (pore volume filled by ice / total pore volume)
526                  saturation_water_ice(ig, ik) = min(ice(ig, ik) / (rho_H2O_ice * porosity_ice_free(ig, ik)), 0.999D-0)
527                 
528                  if (saturation_water_ice(ig, ik).lt.0.D0) then
529                        print *, "WARNING!! soil water ice negative at ", ig, ik
530                        print *, "saturation value: ", saturation_water_ice(ig, ik)
531                        print *, "setting saturation to 0"
532                        saturation_water_ice(ig, ik) = max(saturation_water_ice(ig, ik), 0.D0)
533                  endif
534                 
535                  ! save the old porosity
536                  porosity_prev(ig, ik) = porosity(ig, ik)
537
538                  ! calculate the new porosity
539                  porosity(ig, ik) = porosity_ice_free(ig, ik) * (1.D0 - saturation_water_ice(ig, ik))
540                  ! Note : saturation_water_ice and porosity are computed from the ice content of the previous timestep
541            enddo
542
543            ! calculate the air density in the first subsurface layer
544            rho(ig) = pplev(ig, 1) / (r * ptsoil(ig, 1))
545
546            ! calculate diffusion coefficients at mid- and interlayers (with ice content from the previous timestep)
547            ! Dependence on saturation_water_ice taken from the relation obtained by Hudson et al. (2009) and Meslin et al. (SSSAJ, 2010)
548            do ik = 1, nsoil  ! loop over all soil levels
549
550                  ! calculate H2O mean thermal velocity
551                  vth(ig, ik) = dsqrt(8.D0 * 8.314D0 * ptsoil(ig, ik) / (pi * 18.D-3))
552
553                  ! The diffusion coefficient is calculated
554
555                  ! calculate the H2O - CO2 collision integral (after Mellon and Jakosky, JGR, 1993)
556                  omega(ig, ik) = 1.442D0 - 0.673D0 * dlog(2.516D-3 * ptsoil(ig, ik)) &
557                        + 0.252D0 * (dlog(2.516D-3 * ptsoil(ig, ik))) ** 2.D0 &
558                        - 0.047D0 * (dlog(2.516D-3 * ptsoil(ig, ik))) ** 3.D0 &
559                        + 0.003D0 * (dlog(2.516D-3 * ptsoil(ig, ik))) ** 4.D0
560
561                  ! calculate the molecular diffusion coefficient
562                  Dm(ig, ik) = 4.867D0 * ptsoil(ig, ik) ** (3.D0 / 2.D0) &
563                        / (pplev(ig, 1) * omega(ig, ik)) * 1.D-4
564
565                  ! calculate the Knudsen diffusion coefficient (divided by D0 = porosity / tortuosity in this case)
566                  Dk(ig, ik) = Dk0 / D0(ig, ik) * meshsize(ig, ik) * vth(ig, ik)
567
568                  ! calculate the midlayer diffusion coeff. (with dependence on saturation_water_ice from Meslin et al., 2010 -  only exact for normal diffusion)
569                  D_mid(ig, ik) = D0(ig, ik) * (1.D0 - saturation_water_ice(ig, ik)) ** 2.D0 * 1.D0 / (1.D0 / Dm(ig, ik) + 1.D0 / Dk(ig, ik))
570
571                  ! Midlayer diffusion coeff. (correlation by Rogers and Nielson, 1991)
572                  !                             D_mid(ig, ik) = D0(ig, ik) * (1. - saturation_water_ice(ig, ik)) * exp(-6. * saturation_water_ice(ig, ik)  &
573                  !                                   * porosity_ice_free(ig, ik) - 6. * saturation_water_ice(ig, ik) ** (14. * porosity_ice_free(ig, ik))) &
574                  !                                   * 1. / (1. / Dm(ig, ik) + 1. / Dk(ig, ik))
575            enddo
576
577            ! calculate the interlayer diffusion coefficient as interpolation of the midlayer coefficients
578            do ik = 1, nsoil - 1
579                  Dm_inter(ig, ik) = ( (layer(ik) - mlayer(ik - 1)) * Dm(ig, ik + 1) &
580                        + (mlayer(ik) - layer(ik)) * Dm(ig, ik) ) / (mlayer(ik) - mlayer(ik - 1))
581                       
582                  Dk_inter(ig, ik) = ( (layer(ik) - mlayer(ik - 1)) * Dk(ig, ik + 1) &
583                        + (mlayer(ik) - layer(ik)) * Dk(ig, ik) ) / (mlayer(ik) - mlayer(ik - 1))
584                       
585!                 saturation_water_ice_inter(ig, ik) = ( (layer(ik) - mlayer(ik - 1)) * saturation_water_ice(ig, ik + 1) &
586!                       + (mlayer(ik) - layer(ik)) * saturation_water_ice(ig, ik) ) / (mlayer(ik) - mlayer(ik - 1))
587                 
588                  if (close_bottom(ig, ik).or.close_top(ig, ik+1)) then
589                        saturation_water_ice_inter(ig, ik) = 0.999D0
590                  else
591                        saturation_water_ice_inter(ig, ik) = min(saturation_water_ice(ig, ik + 1), saturation_water_ice(ig, ik))  ! new diffusion interaction
592                  endif
593                 
594                  saturation_water_ice_inter(ig, ik) = min(saturation_water_ice(ig, ik + 1), saturation_water_ice(ig, ik))  ! new diffusion interaction
595                 
596                  D_inter(ig, ik) = D0(ig, ik) * (1.D0 - saturation_water_ice_inter(ig, ik)) ** 2.D0 * 1.D0 / (1.D0 / Dm_inter(ig, ik) + 1.D0 / Dk_inter(ig, ik))
597                 
598!                  D_inter(ig, ik) = ( (layer(ik) - mlayer(ik - 1)) * D_mid(ig, ik + 1) &  ! old implementation with direct interpolation
599!                        + (mlayer(ik) - layer(ik)) * D_mid(ig, ik) ) / (mlayer(ik) - mlayer(ik - 1))
600            enddo
601      endif
602enddo
603
604! Computation of new adsorption / desorption constants (Ka, Kd coefficients), and C, E, F coefficients
605! ===================================================================================================
606
607do ig = 1, ngrid  ! loop over all gridpoints
608      if(.not.watercaptag(ig)) then  ! if there is no polar cap
609            do ik = 1, nsoil  ! loop over all soil levels
610
611                  if (choice_ads.eq.1) then  ! Constant, uniform diffusion coefficient D0 (prescribed)
612
613                        ! first segment of the bilinear function (before monolayer saturation)
614                        ! calculate the equilibrium adsorption coefficient
615                        k_ads_eq(ig, ik) = rho_soil(ig, ik) * S * vth(ig, ik) / 4.D0 &
616                              / (dexp(-enthalpy_ads / (8.314D0 * ptsoil(ig, ik))) / tau0)
617
618                        ! calculate the desorption coefficient
619                        Kd(ig, ik) = kinetic_factor / (1.D0 + k_ads_eq(ig, ik) / porosity(ig, ik))
620
621                        ! calculate the absorption coefficient
622                        Ka(ig, ik) = kinetic_factor * k_ads_eq(ig, ik) / (1.D0 + k_ads_eq(ig, ik) / porosity(ig, ik))
623
624
625                        ! second segment of the bilinear function (after monolayer saturation) ! added 2020
626                        ! calculate the equilibrium adsorption coefficient
627                        k_ads_eq2(ig, ik) = rho_soil(ig, ik) * S * vth(ig, ik) / 4.D0 &
628                              / (dexp(-enthalpy_ads2 / (8.314D0 * ptsoil(ig, ik))) / tau0) ! added 2020
629
630                        ! calculate the desorption coefficient
631                        Kd2(ig, ik) = kinetic_factor / (1.D0 + k_ads_eq2(ig, ik) / porosity(ig, ik)) ! added 2020
632
633                        ! calculate the absorption coefficient
634                        Ka2(ig, ik) = kinetic_factor * k_ads_eq2(ig, ik) / (1.D0 + k_ads_eq2(ig, ik) / porosity(ig, ik)) ! added 2020
635
636                        !                             Kd(ig, ik) = exp(-enthalpy_ads / (8.314 * ptsoil(ig, ik)))
637                        !     &                             / tau0     ! * 1.D-18
638                        !                             Ka(ig, ik) = rho_soil(ig, ik) * S * vth(ig, ik) / 4.   ! * 1.D-18
639
640!                        if ((ngrid.eq.1).and.(ik.eq.18)) then  ! if 1D simulation and uppermost layer
641!                              print * , 'Ka, Kd, Ka / Kd = ', Ka(ig, ik), Kd(ig, ik), Ka(ig, ik) / Kd(ig, ik)
642!                              print * , 'k_ads_eq = ', k_ads_eq(ig, ik)
643!                              print * , 'porosity = ', porosity(ig, ik)
644!                        endif
645
646                  elseif (choice_ads.eq.2) then ! modified 2020 (with DeltaQ instead of Q)
647                        !                                   Kd(ig, ik) = Kd_ref * exp(-enthalpy_ads / 8.314 *
648                        !     &                                   (1. / ptsoil(ig, ik) - 1. / Tref))
649                        !                                   Ka(ig, ik) = rho_soil(ig, ik) * Ka_ref * sqrt(ptsoil(ig, ik) / Tref)
650                        !     &                                   / nref
651
652                        ! first segment of the bilinear function (before monolayer saturation)
653                        ! calculate the equilibrium adsorption coefficient
654                        k_ads_eq(ig, ik) = 8.314D0 * ptsoil(ig,ik)/18.D-3 *Kref * dexp(DeltaQ_ads / 8.314D0 *  &
655                              (1.D0 / ptsoil(ig, ik) - 1.D0 / Tref)) !  Changed enthalpy_ads to DeltaQ_ads to ensure correct dependance/behaviour of k_ads_eq with temperature
656
657                        ! calculate the desorption coefficient
658                        Kd(ig, ik) = kinetic_factor / (1.D0 + k_ads_eq(ig, ik) / porosity(ig, ik))
659
660                        ! calculate the absorption coefficient
661                        Ka(ig, ik) = kinetic_factor * k_ads_eq(ig, ik) / (1.D0 + k_ads_eq(ig, ik) / porosity(ig, ik))
662
663                        ! second segment of the bilinear function (after monolayer saturation) ! added 2020
664                        ! calculate the equilibrium adsorption coefficient
665                        k_ads_eq2(ig, ik) = 8.314D0 * ptsoil(ig,ik)/18.D-3 * Kref2 * dexp(DeltaQ_ads2 / 8.314D0 *  &
666                              (1.D0 / ptsoil(ig, ik) - 1.D0 / Tref)) ! Changed enthalpy_ads2 to DeltaQ_ads2 to ensure correct dependance/behaviour of k_ads_eq with temperature
667
668                        ! calculate the desorption coefficient
669                        Kd2(ig, ik) = kinetic_factor / (1.D0 + k_ads_eq2(ig, ik) / porosity(ig, ik))
670
671                        ! calculate the absorption coefficient
672                        Ka2(ig, ik) = kinetic_factor * k_ads_eq2(ig, ik) / (1.D0 + k_ads_eq2(ig, ik) / porosity(ig, ik))
673                  endif
674                 
675                  ! calculate the amount of water vapor at monolayer saturation ! modified 2020
676                  nsat(ig, ik) = adswater_sat_mono * Kd(ig, ik) / Ka(ig, ik) !
677
678                  ! calculate the intercept of the second line in the bilinear function ! added 2020; corrected 2024
679                  c0(ig, ik) = ( 1.D0 - (k_ads_eq2(ig, ik) / k_ads_eq(ig, ik))) * adswater_sat_mono
680
681                  ! calculate C, E, and F coefficients for later calculations
682                  C(ig, ik) = interlayer_dz(ig, ik) / ptimestep
683
684                  ! first segment of the bilinear function (before monolayer saturation)
685                  E(ig, ik) = Ka(ig, ik) / (1.D0 + ptimestep * Kd(ig, ik))
686                  F(ig, ik) = Kd(ig, ik) / (1.D0 + ptimestep * Kd(ig, ik))
687
688                  ! second segment of the bilinear function (after monolayer saturation) ! added 2020
689                  E2(ig, ik) = Ka2(ig, ik) / (1.D0 + ptimestep * Kd2(ig, ik))
690                  F2(ig, ik) = Kd2(ig, ik) / (1.D0 + ptimestep * Kd2(ig, ik))
691
692                  ! save the old water concentration
693                  znsoilprev(ig, ik) = znsoil(ig, ik)
694                  znsoilprev2(ig, ik) = znsoil(ig, ik)  ! needed in "special case" loop because adswprev can be changed before this loop
695                  adswprev(ig, ik) = adswater(ig, ik)  !!!! Besoin de adswprev2 ???
696                  iceprev(ig, ik) = ice(ig, ik)  ! needed in "special case" loop for the same reason
697                 
698                  ! calculate the total ice + water vapor content
699                  ztot1(ig, ik) = porosity_prev(ig, ik) * znsoil(ig, ik) + ice(ig, ik)
700            enddo
701      endif
702enddo
703
704! Computation of delta, alpha and beta coefficients ETC !!!!
705! ===========================================================
706do ig = 1, ngrid  ! loop over all gridpoints
707      ! initialize the flux to zero to avoid undefined values
708      zdqsdifrego(ig) = 0.D0
709      do ik = 0, nsoil - 1
710            flux(ig, ik) = 0.
711      enddo
712     
713      if(.not.watercaptag(ig)) then  ! if there is no polar cap
714            do ik = 1, nsoil  ! loop over all soil levels
715                  if (ice(ig, ik).eq.0.) then
716                        ice_index_flag(ik) = .false.
717                  else
718                        ice_index_flag(ik) = .true.
719                  endif
720            enddo
721
722
723            do ik = 1, nsoil  ! loop over all soil levels
724                 
725                  ! calculate A and B coefficients ! modified 2020
726                  if ( .not.over_mono_sat_flag(ig, ik) ) then ! Assume cell below the monolayer saturation
727                        ! calculate A and B coefficients (first segment of the bilinear function)
728                        A(ig, ik) = E(ig, ik) * interlayer_dz(ig, ik)
729                        B(ig, ik) = interlayer_dz(ig, ik) * F(ig, ik) * adswprev(ig, ik)
730                  else ! Assume cell over the monolayer saturation ! added 2020
731                        ! calculate A and B coefficients (second segment of the bilinear function)
732                        A(ig, ik) = E2(ig, ik) * interlayer_dz(ig, ik)
733                        B(ig, ik) = interlayer_dz(ig, ik) * F2(ig, ik) * adswprev(ig, ik) &
734                                  + interlayer_dz(ig, ik) * (F2(ig, ik) * Kd2(ig, ik) * c0(ig, ik) * ptimestep - Kd2(ig,ik)*c0(ig,ik)) ! Corrected 2024
735                  endif
736                 
737                  ! calculate the saturation pressure
738                  P_sat_soil(ig, ik) = 611.D0 * dexp(22.5D0 * (1.D0 - 273.16D0 / ptsoil(ig, ik))) ! maybe take a new expression (Hsublim = 51.1 kJ / mol) ?
739                   
740                  ! calculate the gas number density at saturation pressure
741                  nsatsoil(ig, ik) = (P_sat_soil(ig, ik) * mw) / (kb * ptsoil(ig, ik))
742            enddo
743
744            ! calculate the alpha, beta, and delta coefficients for the surface layer
745            delta0(ig) = pa(ig, 1) + pb(ig, 2) * (1.D0 - pd(ig, 2)) + porosity_ice_free(ig, 1) * pb(ig, 1)
746            alpha0(ig) = porosity_ice_free(ig, 1) * pb(ig, 1) / delta0(ig)
747           
748            beta0(ig) = ( pb(ig, 2) * pc(ig, 2) + pq(ig, 1, igcm_h2o_vap) * pa(ig, 1) + pqsurf(ig) ) &
749                  / delta0(ig)
750                 
751            ! set loop flag
752            do ik = 1, nsoil
753                  condensation_flag = .false.
754            enddo
755           
756            sublimation_flag = .true.
757            sublim_n = 0
758            do while (sublimation_flag)  ! while there is sublimation
759!                  print *, "sublimation loop: ", sublim_n
760                  sublim_n = sublim_n + 1
761                 
762                  if (sublim_n.gt.100) then
763                        print *, "sublimation not converging in call ", n
764                        print *, "might as well stop now"
765                        stop
766                  endif
767                 
768                  ! reset flag
769                  sublimation_flag = .false.
770                 
771                  recompute_all_cells_ads_flag = .true.
772                  satur_mono_n = 0
773                  do while (recompute_all_cells_ads_flag)
774!                        print *, satur_mono_n
775                        satur_mono_n = satur_mono_n + 1
776                       
777                        ! reset loop flag
778                        recompute_all_cells_ads_flag = .false.
779                       
780                        ! calculate the surface air density
781                        rhosurf(ig) = pplev(ig, 1) / (r * ptsrf(ig))
782                       
783                        ! calculate the coefficients in the upermost layer
784                        if (exchange(ig)) then  ! if there is surface exchange
785                             
786                              ! calculate the delta, alpha, and beta coefficients
787                              zdelta(ig, 1) = A(ig, 1) + porosity(ig, 1) * C(ig, 1) &
788                                    + D_inter(ig, 1) / midlayer_dz(ig, 1) &  !!! est - ce que le terme pb inclut le bon rho ?
789                                    + porosity_ice_free(ig, 1) * pb(ig, 1) / (rho(ig) * ptimestep) * (1.D0 - alpha0(ig))
790                                    !pourrait remplacer pb/(ptime*rho(1/2)) par zcdv/ptime
791                             
792                              zalpha(ig, 1) = D_inter(ig, 1) / midlayer_dz(ig, 1) * 1.D0 / zdelta(ig, 1)
793                             
794                              zbeta(ig, 1) = ( porosity_ice_free(ig, 1) * pb(ig, 1) / ptimestep * beta0(ig) &
795                                    + porosity_prev(ig, 1) * C(ig, 1) * znsoilprev(ig, 1) + B(ig, 1) ) &
796                                    / zdelta(ig, 1)
797                        else
798                              ! calculate the delta, alpha, and beta coefficients without exchange
799                              zdelta(ig, 1) = A(ig, 1) + D_inter(ig, 1) / midlayer_dz(ig, 1)  &
800                                    + D_mid(ig, 1) / midlayer_dz(ig, 0) &
801                                    + porosity(ig, 1) * C(ig, 1)
802                                   
803                              zalpha(ig, 1) = D_inter(ig, 1) / midlayer_dz(ig, 1) * 1. / zdelta(ig, 1)
804                             
805                              zbeta(ig, 1) = ( D_mid(ig, 1) / midlayer_dz(ig, 0) * qsat_surf(ig) * rhosurf(ig) &
806                                    + porosity_prev(ig, 1) * C(ig, 1) * znsoilprev(ig, 1) + B(ig, 1) ) &
807                                    / zdelta(ig, 1)
808                                       
809                        endif
810                       
811                        ! check for ice   
812                        if (ice_index_flag(1)) then
813                              ! set the alpha coefficient to zero
814                              zalpha(ig, 1) = 0.D0
815                              ! set the beta coefficient to the number density at saturation pressure
816                              zbeta(ig, 1) = nsatsoil(ig, 1)
817                        endif
818
819                        do ik = 2, nsoil - 1  ! go through all other layers
820                       
821                              ! calculate delta, alpha, and beta coefficients
822                             
823                              zdelta(ig, ik) = A(ig, ik) + porosity(ig, ik) * C(ig, ik) &
824                                    + D_inter(ig, ik) / midlayer_dz(ig, ik) &
825                                    + D_inter(ig, ik - 1) / midlayer_dz(ig, ik - 1) * (1.D0 - zalpha(ig, ik - 1))
826
827                              zalpha(ig, ik) = D_inter(ig, ik) / midlayer_dz(ig, ik) * 1.D0 / zdelta(ig, ik)
828
829                              zbeta(ig, ik) = ( D_inter(ig, ik - 1) / midlayer_dz(ig, ik - 1) * zbeta(ig, ik - 1) &
830                                    + B(ig, ik) &
831                                    + porosity_prev(ig, ik) * C(ig, ik) * znsoilprev(ig, ik) ) / zdelta(ig, ik)
832                             
833                              ! check for ice   
834                              if (ice_index_flag(ik)) then
835                                    ! set the alpha coefficient to zero
836                                    zalpha(ig, ik) = 0.D0
837                                    ! set the beta coefficient to the number density at saturation pressure
838                                    zbeta(ig, ik) = nsatsoil(ig, ik)
839                              endif
840                        enddo
841
842                        ! Computation of pqsoil, ztot1, zq1temp2 and zdqsdifrego
843                        ! =======================================================  =
844                       
845                        ! calculate the preliminary amount of water vapor in the bottom layer
846                        if (ice_index_flag(nsoil)) then  ! if there is ice
847                              ! set the vapor number density to saturation
848                              znsoil(ig, nsoil) = nsatsoil(ig, nsoil)
849                        else
850                              ! calculate the vapor number density in the lowest layer
851                              znsoil(ig, nsoil) = ( D_inter(ig, nsoil - 1) / midlayer_dz(ig, nsoil - 1) * zbeta(ig, nsoil - 1) &
852                                    + porosity_prev(ig, nsoil) * C(ig, nsoil) * znsoilprev(ig, nsoil) + B(ig, nsoil) ) &
853                                    / ( porosity(ig, nsoil) * C(ig, nsoil) + A(ig, nsoil) &
854                                    + D_inter(ig, nsoil - 1) / midlayer_dz(ig, nsoil - 1) * (1.D0 - zalpha(ig, nsoil - 1)) )
855                        endif
856                       
857                        ! calculate the preliminary amount of adsorbed water
858                        if (.not.over_mono_sat_flag(ig, nsoil)) then ! modified 2024
859                            adswater_temp(ig, nsoil) = ( Ka(ig, nsoil) * ptimestep * znsoil(ig, nsoil) + adswprev(ig, nsoil) ) &
860                                 / (1.D0 + ptimestep * Kd(ig, nsoil))
861                        else
862                             adswater_temp(ig, nsoil) = ( Ka2(ig, nsoil) * ptimestep * znsoil(ig, nsoil) + adswprev(ig, nsoil) + ptimestep * c0(ig,nsoil) * Kd2(ig,nsoil)) &
863                                 / (1.D0 + ptimestep * Kd2(ig, nsoil))
864                        endif                             
865
866
867 
868                        do ik = nsoil-1, 1, -1  ! backward loop over all layers
869                                 
870                              znsoil(ig, ik) = zalpha(ig, ik) * znsoil(ig, ik + 1) + zbeta(ig, ik)
871                             
872                              if (.not.over_mono_sat_flag(ig, nsoil)) then ! modified 2024
873                                adswater_temp(ig, ik) = ( Ka(ig, ik) * ptimestep * znsoil(ig, ik) + adswprev(ig, ik) ) &
874                                    / (1.D0 + ptimestep * Kd(ig, ik))
875                              else
876                                 adswater_temp(ig, ik) = ( Ka2(ig, ik) * ptimestep * znsoil(ig, ik) + adswprev(ig, ik) + ptimestep * c0(ig,ik) * Kd2(ig,ik)) &
877                                    / (1.D0 + ptimestep * Kd2(ig, ik))
878                              endif 
879                        enddo
880                       
881                        ! check if any cell is over monolayer saturation
882                        do ik = 1, nsoil  ! loop over all soil levels
883
884                              ! if( (ik.lt.nsoil) .and. (recompute_cell_ads_flag(ik+1) = .true.) ) then ! Make this loop faster as all cells also need to be recomputed if the cell below needs to be recomputed ! ARNAU
885                              !     recompute_cell_ads_flag(ik) = .true.
886                              !     cycle ! Jump loop
887                              ! endif
888                             
889                              ! Change of coefficients ! ARNAU
890                              recompute_cell_ads_flag(ik) = .false. ! Assume there is no change of coefficients
891                              if ( adswater_temp(ig, ik).ge.adswater_sat_mono ) then
892
893                                    print *, "NOTICE: over monolayer saturation"
894
895                                    if ( .not.over_mono_sat_flag(ig, ik) ) then
896                                          ! If the cell enters this scope, it
897                                          ! means that the cell is over the monolayer
898                                          ! saturation after calculation, but the
899                                          ! coefficients assume it is below. Therefore,
900                                          ! the cell needs to be recomputed
901                                   
902                                          over_mono_sat_flag(ig, ik) = .true.
903                                          recompute_cell_ads_flag(ik) = .true. ! There is a change of coefficients
904
905                                          ! change the A and B coefficients (change from first segment to second segment of the bilinear function, as we are over the monolayer saturation is_cell_over_monolayer_sat)
906                                          A(ig, ik) = E2(ig, ik) * interlayer_dz(ig, ik)
907                                          B(ig, ik) = interlayer_dz(ig, ik) * F2(ig, ik) * adswprev(ig, ik) &
908                                                          + interlayer_dz(ig, ik) * (F2(ig, ik) * Kd2(ig, ik) * c0(ig, ik) * ptimestep - Kd2(ig,ik)*c0(ig,ik)) ! Corrected 2024
909                                    endif
910                              else
911                                    if ( over_mono_sat_flag(ig, ik) ) then
912                                          ! If the cell enters this scope, it
913                                          ! means that the cell is below the monolayer
914                                          ! saturation after calculation, but the
915                                          ! coefficients assume it is above. Therefore,
916                                          ! the cell needs to be recomputed
917                                   
918                                          over_mono_sat_flag(ig, ik) = .false.
919                                          recompute_cell_ads_flag(ik) = .true. ! There is a change of coefficients
920
921                                          ! calculate A and B coefficients (reset to first segment of the bilinear function)
922                                          A(ig, ik) = E(ig, ik) * interlayer_dz(ig, ik)
923                                          B(ig, ik) = interlayer_dz(ig, ik) * F(ig, ik) * adswprev(ig, ik)
924                                    endif
925                              endif
926                        enddo
927                       
928                        ! if one cell needs to be recomputed, then all the column is to be recomputed too
929                        do ik = 1, nsoil
930                              if ( recompute_cell_ads_flag(ik) ) then
931                                    recompute_all_cells_ads_flag = .true.
932                              else
933                                    adswater(ig, ik) = adswater_temp(ig, ik) ! if no recomputing is needed, store the value (it may be wrong if the cell below needs to be recomputed. It will be correct in the next loop iterations)
934                              endif
935                        enddo
936                  enddo  ! end loop while recompute_all_cells_ads_flag (monolayer saturation)
937                 
938                  !? I'm not sure if this should be calculated here again. I have a feeling that ztot1 is meant
939                  ! as the value from the previous timestep
940                  !do ik = 1, nsoil
941                  !      ztot1(ig, ik) = porosity_prev(ig, ik) * znsoil(ig, ik) + ice(ig, ik)
942                  !enddo
943
944                  ! Calculation of Fnet + Fads
945                  ! =============================
946                 
947                  ! calulate the flux variable
948                 
949                  ! calculate the flux in the top layer
950                  if (exchange(ig)) then  ! if there is exchange
951                        ! calculate the flux taking diffusion to the atmosphere into account
952                        flux(ig, 0) = porosity_ice_free(ig, 1) * pb(ig, 1) / ptimestep &
953                              * ( znsoil(ig, 1) / rho(ig) - beta0(ig) - alpha0(ig) * znsoil(ig, 1) / rho(ig) )   
954                  elseif (pqsurf(ig).gt.0.) then
955                        ! assume that the surface is covered by waterice (if it is co2ice it should not call this subroutine at all)
956                        flux(ig, 0) = D_mid(ig, 1) / midlayer_dz(ig, 0) * ( znsoil(ig, 1) - qsat_surf(ig) * rhosurf(ig) )
957                  endif
958                 
959                  ! check if the ground would take up water from the surface but there is none
960                  if ((.not.exchange(ig)).and.(pqsurf(ig).eq.0.).and.(flux(ig, 0).lt.0.)) then
961                        flux(ig, 0) = 0.
962                  endif
963                 
964                  ! calculate the flux in all other layers
965                  do ik = 1, nsoil - 1
966                        flux(ig, ik) = D_inter(ig, ik) * ( znsoil(ig, ik + 1) - znsoil(ig, ik) ) / midlayer_dz(ig, ik)
967                  enddo
968                 
969                  ! calculate dztot1 which descibes the change in ztot1 (water vapor and ice). It is only used for output (directly and indirectly)
970                 
971                  ! calculate the change in ztot1             
972                 
973                  do ik = 1, nsoil - 1
974                        dztot1(ik) = ( flux(ig, ik) - flux(ig, ik - 1) ) / interlayer_dz(ig, ik) &
975                              - A(ig, ik) / interlayer_dz(ig, ik) * znsoil(ig, ik) &
976                              + B(ig, ik) / interlayer_dz(ig, ik)
977                  enddo
978                 
979                  dztot1(nsoil) = -flux(ig, nsoil - 1) / interlayer_dz(ig, nsoil) &
980                        - A(ig, nsoil) / interlayer_dz(ig, nsoil) * znsoil(ig, nsoil) &
981                        + B(ig, nsoil) / interlayer_dz(ig, nsoil)
982                 
983                  ! Condensation / sublimation
984                  do ik = 1, nsoil  ! loop over all layers
985!                        print *, "ice in layer      ", ik, ": ", ice(ig, ik)
986!                        print *, "vapor in layer    ", ik, ": ", znsoil(ig, ik)
987!                        print *, "sat dens in layer ", ik, ": ", nsatsoil(ig, ik)
988!                        print *, ""
989                        if (ice_index_flag(ik)) then  ! if there is ice
990                             
991                              ! Compute ice content
992                              ice(ig, ik) = ztot1(ig, ik) + dztot1(ik) * ptimestep - porosity(ig, ik) * nsatsoil(ig, ik)
993
994                              if (ice(ig, ik).lt.0.D0) then  ! if all the ice is used up
995!                                    print *, "########complete sublimation in layer", ik, "  cell:", ig
996!                                    print *, "ztot1: ", ztot1(ig, ik)
997!                                    print *, "dztot1*timestep: ", dztot1(ik) * ptimestep
998!                                    print *, "vapour: ", porosity(ig, ik) * nsatsoil(ig, ik)
999!                                    print *, "znsoil: ", znsoil(ig, ik)
1000!                                    print *, "nsatsoil: ", nsatsoil(ig, ik)
1001!                                    print *, "porosity: ", porosity(ig, ik)
1002!                                    print *, "ice: ", ice(ig, ik)
1003!                                    print *, "exchange: ", exchange(ig)
1004!                                    print *, "co2ice: ", co2ice(ig)
1005!                                    print *, ""
1006!                                    print *, "zalpha: ", zalpha(ig, ik)
1007!                                    print *, "zbeta: ", zbeta(ig, ik)
1008!                                    print *, ""
1009                                   
1010                                    ! set the ice content to zero
1011                                    ice(ig, ik) = 0.D0
1012                                   
1013                                    ! change the water vapor from the previous timestep. (Watch out! could go wrong)
1014                                    znsoilprev(ig, ik) = ztot1(ig, ik) / porosity_prev(ig, ik)
1015!                                    print *, "new znsoilprev", znsoilprev(ig, ik)
1016                                   
1017                                    ! remove the ice flag and raise the sublimation flag
1018                                    ice_index_flag(ik) = .false.
1019                                    sublimation_flag = .true.
1020                              endif
1021
1022                        elseif (znsoil(ig, ik).gt.nsatsoil(ig, ik)) then  ! if there is condenstation
1023                              !ice(ig, ik) = ztot1(ig, ik) + dztot1(ik) * ptimestep - porosity(ig, ik) * nsatsoil(ig, ik)
1024
1025                             
1026!                              print *, "=========== new condensation in layer", ik, "  cell:", ig
1027!                              print *, "znsoil, nsatsoil: ", znsoil(ig, ik), nsatsoil(ig, ik)
1028!                              print *, "ztot1: ", ztot1(ig, ik)
1029!                              print *, "dztot1: ", dztot1(ik)
1030!                              print *, "ice: ", ice(ig,ik)
1031!                              print *, ""
1032!                              print *, "zalpha: ", zalpha(ig, ik)
1033!                              print *, "zbeta: ", zbeta(ig, ik)
1034!                              print *, ""
1035                             
1036                              !if (ice(ig, ik).lt.0.D0) then
1037                                    ! set the ice content to zero
1038                              !      ice(ig, ik) = 0.D0
1039                              if (.not.condensation_flag(ik)) then
1040                                    ! raise the ice and sublimation flags
1041                                    ice_index_flag(ik) = .true.
1042                                    sublimation_flag = .true.
1043!                                    print *, condensation_flag(ik)
1044                                    condensation_flag(ik) = .true.
1045!                                    print *, condensation_flag(ik)
1046!                                    print *, "set condensation flag"
1047!                              else
1048!                                    print *, "condensation loop supressed"
1049                              endif
1050                        endif
1051                       
1052                        ! calculate the difference between total water content and ice + vapor content (only used for output)
1053                        diff(ig, ik) = ice(ig, ik) + porosity(ig, ik) * znsoil(ig, ik) &
1054                              + adswater(ig, ik) - ztot1(ig, ik) - dztot1(ik) * ptimestep
1055                  enddo  ! loop over all layer
1056            enddo  ! end first loop while sublimation_flag (condensation / sublimation)
1057
1058            if (exchange(ig)) then  ! if there is exchange
1059                 
1060                  ! calculate the temporty mixing ratio above the surface
1061                  zq1temp2(ig) = beta0(ig) + alpha0(ig) * znsoil(ig, 1) / rho(ig)
1062                  ! calculate the flux from the subsurface
1063                  zdqsdifrego(ig) = porosity_ice_free(ig, 1) * pb(ig, 1) / ptimestep * (zq1temp2(ig) - znsoil(ig, 1) / rho(ig) )
1064                 
1065            else
1066                  ! calculate the flux from the subsurface
1067                  zdqsdifrego(ig) = D_mid(ig, 1) / midlayer_dz(ig, 0) * (znsoil(ig, 1) - qsat_surf(ig) * rhosurf(ig)) ! faut - il faire intervenir la porosite ?
1068                  if ((pqsurf(ig).eq.0.).and.(zdqsdifrego(ig).lt.0.)) then
1069                        zdqsdifrego(ig) = 0.
1070                  endif
1071            endif
1072
1073! Special case where there is not enough ice on the surface to supply the subsurface for the whole timestep
1074! (when exchange with the atmosphere is disabled): the upper boundary condition becomes a flux condition
1075! (and not qsat_surf) and all the remaining surface ice is sublimated and transferred to the subsurface.
1076! the actual change in surface ice happens in a higher subroutine
1077! =========================================================================================================  =
1078
1079            if ( (.not.exchange(ig)) &
1080            .and. ( (-zdqsdifrego(ig) * ptimestep) &
1081            .gt.( pqsurf(ig) + pdqsdifpot(ig, igcm_h2o_ice) * ptimestep) ) &
1082            .and.( (pqsurf(ig) + pdqsdifpot(ig, igcm_h2o_ice) * ptimestep).gt.0. ) ) then
1083
1084                  ! calculate a new flux from the subsurface
1085                  zdqsdifrego(ig) = -( pqsurf(ig) + pdqsdifpot(ig, igcm_h2o_ice) * ptimestep ) / ptimestep
1086                 
1087!                  ! check case where there is CO2 ice on the surface but qsurf = 0
1088!                  if (co2ice(ig).gt.0.) then
1089!                        zdqsdifrego(ig) = 0.D0 
1090!                  endif
1091
1092
1093                  ! calculate A and B coefficients ! modified 2020
1094                  if ( .not.over_mono_sat_flag(ig, ik) ) then ! Assume cell below the monolayer saturation
1095                        ! calculate A and B coefficients (first segment of the bilinear function)
1096                        A(ig, ik) = E(ig, ik) * interlayer_dz(ig, ik)
1097                        B(ig, ik) = interlayer_dz(ig, ik) * F(ig, ik) * adswprev(ig, ik)
1098                  else ! Assume cell over the monolayer saturation ! added 2020
1099                        ! calculate A and B coefficients (second segment of the bilinear function)
1100                        A(ig, ik) = E2(ig, ik) * interlayer_dz(ig, ik)
1101                        B(ig, ik) = interlayer_dz(ig, ik) * F2(ig, ik) * adswprev(ig, ik) &
1102                                  + interlayer_dz(ig, ik) * (F2(ig, ik) * Kd2(ig, ik) * c0(ig, ik) * ptimestep - Kd2(ig,ik)*c0(ig,ik)) ! Corrected 2024
1103                  endif
1104                 
1105                  ! initialize all flags for the loop
1106                  do ik = 1, nsoil
1107                       
1108                        ! initialize the ice
1109                        ice(ig, ik) = iceprev(ig, ik)
1110                        if (iceprev(ig, ik).eq.0.) then
1111                              ice_index_flag(ik) = .false.
1112                        else
1113                              ice_index_flag(ik) = .true.
1114                        endif
1115                  enddo
1116
1117                  ! loop while there is new sublimation
1118                 
1119                  sublimation_flag = .true.
1120                  sublim_n = 0
1121                  do while (sublimation_flag)
1122!                        print *, "special case sublimation loop: ", sublim_n
1123                        sublim_n = sublim_n + 1
1124                        if (sublim_n.gt.100) then
1125                              print *, "special case sublimation not converging in call ", n
1126                              print *, "might as well stop now"
1127                              stop
1128                        endif
1129                       
1130                        ! reset the sublimation flag
1131                        sublimation_flag = .false.
1132                                               
1133                        ! loop until good values accounting for monolayer saturation are found
1134                        recompute_all_cells_ads_flag = .true.
1135                        do while (recompute_all_cells_ads_flag)
1136                             
1137                              ! reset loop flag
1138                              recompute_all_cells_ads_flag = .false.
1139                             
1140                              ! calculate the Delta, Alpha, and Beta coefficients in the top layer
1141                              zdelta(ig, 1) = A(ig, 1) + porosity(ig, 1) * C(ig, 1) + D_inter(ig, 1) / midlayer_dz(ig, 1)
1142                             
1143                              zalpha(ig, 1) = D_inter(ig, 1) / midlayer_dz(ig, 1) * 1.D0 / zdelta(ig, 1)
1144                             
1145                              zbeta(ig, 1) = ( porosity_prev(ig, 1) * C(ig, 1) * znsoilprev2(ig, 1) + B(ig, 1) - zdqsdifrego(ig) ) &
1146                                    / zdelta(ig, 1)
1147                                   
1148                              ! check for ice
1149                              if (ice_index_flag(1)) then
1150                                          ! set the alpha coefficient to zero
1151                                          zalpha(ig, 1) = 0.D0
1152                                          ! set the beta coefficient to the number density at saturation pressure
1153                                          zbeta(ig, 1) = nsatsoil(ig, 1)
1154                              endif
1155
1156                              do ik = 2, nsoil - 1  ! loop over all other layers
1157                             
1158                                    ! calculate the Delta, Alpha, and Beta coefficients in the layer
1159                                    zdelta(ig, ik) = A(ig, ik) + porosity(ig, ik) * C(ig, ik) + D_inter(ig, ik) / midlayer_dz(ig, ik) &
1160                                          + D_inter(ig, ik - 1) / midlayer_dz(ig, ik - 1) * (1.D0 - zalpha(ig, ik - 1))
1161                                                                             
1162                                    zalpha(ig, ik) = D_inter(ig, ik) / midlayer_dz(ig, ik) * 1.D0 / zdelta(ig, ik)
1163                                   
1164                                    zbeta(ig, ik) = ( D_inter(ig, ik - 1) / midlayer_dz(ig, ik - 1) * zbeta(ig, ik - 1) &
1165                                          + porosity_prev(ig, ik) * C(ig, ik) * znsoilprev2(ig, ik) + B(ig, ik) ) / zdelta(ig, ik)       
1166                                   
1167                                    ! check for ice
1168                                    if (ice_index_flag(ik)) then
1169                                          ! set the alpha coefficient to zero
1170                                          zalpha(ig, ik) = 0.D0
1171                                          ! set the beta coefficient to the number density at saturation pressure
1172                                          zbeta(ig, ik) = nsatsoil(ig, ik)
1173                                    endif
1174                              enddo
1175                             
1176                              ! calculate the preliminary amount of water vapor in the bottom layer
1177                              if (ice_index_flag(nsoil)) then
1178                                    ! set the vapor number density to saturation
1179                                    znsoil(ig, nsoil) = nsatsoil(ig, nsoil) 
1180                              else
1181                                    ! calculate the vapor number density in the lowest layer
1182                                    znsoil(ig, nsoil) = ( D_inter(ig, nsoil - 1) / midlayer_dz(ig, nsoil - 1) * zbeta(ig, nsoil - 1) &
1183                                          + porosity_prev(ig, nsoil) * C(ig, nsoil) * znsoilprev2(ig, nsoil) + B(ig, nsoil) ) &
1184                                          / ( porosity(ig, nsoil) * C(ig, nsoil) &
1185                                          + D_inter(ig, nsoil - 1) / midlayer_dz(ig, nsoil - 1) * (1.D0 - zalpha(ig, nsoil - 1)) &
1186                                          + A(ig, nsoil) )
1187                              endif
1188                             
1189
1190                             ! calculate the preliminary amount of adsorbed water
1191                            if (.not.over_mono_sat_flag(ig, nsoil)) then ! modified 2024
1192                                adswater_temp(ig, nsoil) = ( Ka(ig, nsoil) * ptimestep * znsoil(ig, nsoil) + adswprev(ig, nsoil) ) &
1193                                     / (1.D0 + ptimestep * Kd(ig, nsoil))
1194                            else
1195                                 adswater_temp(ig, nsoil) = ( Ka2(ig, nsoil) * ptimestep * znsoil(ig, nsoil) + adswprev(ig, nsoil) + ptimestep * c0(ig,nsoil) * Kd2(ig,nsoil)) &
1196                                 / (1.D0 + ptimestep * Kd2(ig, nsoil))
1197                            endif
1198
1199
1200
1201                            do ik = nsoil-1, 1, -1  ! backward loop over all layers
1202   
1203                                  znsoil(ig, ik) = zalpha(ig, ik) * znsoil(ig, ik + 1) + zbeta(ig, ik)
1204   
1205                                  if (.not.over_mono_sat_flag(ig, nsoil)) then ! modified 2024
1206                                    adswater_temp(ig, ik) = ( Ka(ig, ik) * ptimestep * znsoil(ig, ik) + adswprev(ig, ik) ) &
1207                                        / (1.D0 + ptimestep * Kd(ig, ik))
1208                                  else
1209                                     adswater_temp(ig, ik) = ( Ka2(ig, ik) * ptimestep * znsoil(ig, ik) + adswprev(ig, ik) + ptimestep * c0(ig,ik) * Kd2(ig,ik)) &
1210                                        / (1.D0 + ptimestep * Kd2(ig, ik))
1211                                  endif
1212                            enddo
1213
1214                             
1215                              ! check for any saturation
1216                              do ik = 1, nsoil
1217                                    ! Change of coefficients ! ARNAU
1218                                    recompute_cell_ads_flag(ik) = .false. ! Assume there is no change of coefficients ! ARNAU
1219                                    if ( adswater_temp(ig, ik).gt.adswater_sat_mono ) then
1220                                   
1221                                          print *, "NOTICE: over monolayer saturation"
1222
1223                                          if ( .not.over_mono_sat_flag(ig, ik) ) then
1224                                               ! If the cell enters this scope, it
1225                                               ! means that the cell is over the monolayer
1226                                               ! saturation after calculation, but the
1227                                               ! coefficients assume it is below. Therefore,
1228                                               ! the cell needs to be recomputed
1229                                         
1230                                                over_mono_sat_flag(ig, ik) = .true.
1231                                                recompute_cell_ads_flag(ik) = .true. ! There is a change of coefficients
1232
1233                                                ! change the A and B coefficients (change from first segment to second segment of the bilinear function, as we are over the monolayer saturation is_cell_over_monolayer_sat)
1234                                                A(ig, ik) = E2(ig, ik) * interlayer_dz(ig, ik)
1235                                                B(ig, ik) = interlayer_dz(ig, ik) * F2(ig, ik) * adswprev(ig, ik) &
1236                                                          + interlayer_dz(ig, ik) * (F2(ig, ik) * Kd2(ig, ik) * c0(ig, ik) * ptimestep - Kd2(ig,ik)*c0(ig,ik)) ! Corrected 2024
1237                                          endif
1238                                    else
1239                                          if ( over_mono_sat_flag(ig, ik) ) then
1240                                               ! If the cell enters this scope, it
1241                                               ! means that the cell is below the monolayer
1242                                               ! saturation after calculation, but the
1243                                               ! coefficients assume it is above. Therefore,
1244                                               ! the cell needs to be recomputed
1245                                         
1246                                                over_mono_sat_flag(ig, ik) = .false.
1247                                                recompute_cell_ads_flag(ik) = .true. ! There is a change of coefficients
1248
1249                                                ! calculate A and B coefficients (reset to first segment of the bilinear function)
1250                                                A(ig, ik) = E(ig, ik) * interlayer_dz(ig, ik)
1251                                                B(ig, ik) = interlayer_dz(ig, ik) * F(ig, ik) * adswprev(ig, ik)
1252                                          endif
1253                                    endif
1254                              enddo
1255                             
1256                              ! raise the saturation flag if any layer has saturated and reset the first layer saturation flag
1257                              do ik = 1, nsoil ! modified 2020
1258                                    if ( recompute_cell_ads_flag(ik) ) then
1259                                          recompute_all_cells_ads_flag = .true.
1260                                    else
1261                                          adswater(ig, ik) = adswater_temp(ig, ik) ! if no recomputing is needed, store the value (it may be wrong if the cell below needs to be recomputed. It will be correct in the next loop iterations)
1262                                    endif
1263                              enddo
1264                        enddo  ! end while loop (adsorption saturation)
1265
1266                        ! set the flux to the previously calculated value for the top layer
1267                        flux(ig, 0) = zdqsdifrego(ig)
1268                       
1269                        ! calculate the flux for all other layers
1270                        do ik = 1, nsoil - 1
1271                              flux(ig, ik) = D_inter(ig, ik) * (znsoil(ig, ik + 1) - znsoil(ig, ik)) / midlayer_dz(ig, ik)
1272                        enddo
1273                       
1274                        ! calculate the change in ztot1
1275                        do ik = 1, nsoil - 1
1276                              dztot1(ik) = (flux(ig, ik) - flux(ig, ik - 1)) / interlayer_dz(ig, ik) &
1277                                    - A(ig, ik) / interlayer_dz(ig, ik) * znsoil(ig, ik) &
1278                                    + B(ig, ik) / interlayer_dz(ig, ik)
1279                        enddo
1280
1281                        dztot1(nsoil) = -flux(ig, nsoil - 1) / interlayer_dz(ig, nsoil) &
1282                              - A(ig, nsoil) / interlayer_dz(ig, nsoil) * znsoil(ig, nsoil) &
1283                              + B(ig, nsoil) / interlayer_dz(ig, nsoil)
1284
1285
1286                        ! Condensation / sublimation
1287                        do ik = 1, nsoil
1288                              if (ice_index_flag(ik)) then
1289
1290                                    ! Compute ice content
1291                                    ice(ig, ik) = ztot1(ig, ik) + dztot1(ik) * ptimestep &
1292                                          - porosity(ig, ik) * nsatsoil(ig, ik)
1293
1294                                    if (ice(ig, ik).lt.0.D0) then  ! if all the ice is used up
1295                                   
1296!                                          print *, "############complete sublimation in layer ", ik
1297!                                          print *, "ztot1: ", ztot1(ig, ik)
1298!                                          print *, "dztot1*timestep: ", dztot1(ik) * ptimestep
1299!                                          print *, "vapour: ", porosity(ig, ik) * nsatsoil(ig, ik)
1300!                                          print *, "znsoil: ", znsoil(ig, ik)
1301!                                          print *, "nsatsoil: ", nsatsoil(ig, ik)
1302!                                          print *, "porosity: ", porosity(ig, ik)
1303!                                          print *, "ice: ", ice(ig, ik)
1304!                                          print *, "exchange: ", exchange(ig)
1305!                                          print *, "co2ice: ", co2ice(ig)
1306!                                          print *, ""
1307!                                          print *, "zalpha: ", zalpha(ig, ik)
1308!                                          print *, "zbeta: ", zbeta(ig, ik)
1309!                                          print *, ""
1310                                         
1311                                          ! set the ice content to zero
1312                                          ice(ig, ik) = 0.D0
1313                                         
1314                                         
1315                                          if (znsoil(ig, ik).gt.nsatsoil(ig, ik)) then
1316                                                print *, "WARNING! complete sublimation, but vapor is oversaturated"
1317                                                print *, "special case loop in cell", ig, ik ,"will be supressed"
1318                                          else
1319                                                ! change the water vapor from the previous timestep. (Watch out! could go wrong)
1320                                                znsoilprev2(ig, ik) = ztot1(ig, ik) / porosity_prev(ig, ik)
1321                                                ! print *, "new znsoilprev", znsoilprev(ig, ik)
1322                                               
1323                                                ! remove the ice flag and raise the sublimation flag
1324                                                ice_index_flag(ik) = .false.
1325                                                sublimation_flag = .true.
1326                                          endif
1327                                    endif
1328
1329                              elseif (znsoil(ig, ik).gt.nsatsoil(ig, ik)) then  ! if there is new ice through condensation
1330                                    if (.not.condensation_flag(ik)) then
1331                                    ! raise the ice and sublimation flags
1332                                    ice_index_flag(ik) = .true.
1333                                    sublimation_flag = .true.
1334!                                    print *, condensation_flag(ik)
1335                                    condensation_flag(ik) = .true.
1336!                                    print *, condensation_flag(ik)
1337!                                    print *, "set condensation flag"
1338!                              else
1339!                                    print *, "special case condensation loop supressed"
1340                              endif
1341                              endif
1342                        enddo
1343                  enddo  ! end first loop while (condensation / sublimation)
1344            endif  ! Special case for all ice on the surface is used up
1345
1346            do ik = 1, nsoil
1347
1348                  diff(ig, ik) = ice(ig, ik) + porosity(ig, ik) * znsoil(ig, ik) &  ! only used for output
1349                        + adswater(ig, ik) - adswprev(ig, ik) &
1350                        - ztot1(ig, ik) - dztot1(ik) * ptimestep  !? This is inconsistent, in the not special case the previous adsorbed water is not counted
1351                                                                  !? also this just overwrites the previous calculation if I see that correctly
1352                 
1353                  ! calculate the total amount of water
1354                  ztot1(ig, ik) = porosity(ig, ik) * znsoil(ig, ik) + ice(ig, ik)
1355                  h2otot(ig, ik) = adswater(ig, ik) + ztot1(ig, ik)
1356            enddo
1357      endif  ! if there is no polar cap
1358enddo  ! loop over all gridpoints
1359
1360! check for choking condition
1361do ig = 1, ngrid
1362      if(.not.watercaptag(ig)) then
1363            do ik = 1, nsoil - 1
1364                  if (ice(ig, ik) / (rho_H2O_ice * porosity_ice_free(ig, ik)).gt.choke_fraction) then  ! if the ice is over saturation or choke_fraction
1365                        if ( flux(ig, ik - 1).gt.0.) then
1366                              if (.not.close_top(ig, ik).and.print_closure_warnings) then
1367                                    print *, "NOTICE: closing top of soil layer due to inward flux", ig, ik     
1368                              endif
1369                              close_top(ig, ik) = .true.
1370                        elseif (flux(ig, ik - 1).lt.0.) then
1371                              if (close_top(ig, ik).and.print_closure_warnings) then
1372                                    print *, "NOTICE: opening top of soil layer due to outward flux", ig, ik     
1373                              endif
1374                              close_top(ig, ik) = .false.
1375                        endif
1376                       
1377                        if ( flux(ig, ik).lt.0.) then
1378                              if (.not.close_bottom(ig, ik).and.print_closure_warnings) then
1379                                    print *, "NOTICE: closing bottom of soil layer due to inward flux", ig, ik     
1380                              endif
1381                              close_bottom(ig, ik) = .true.
1382                        elseif (flux(ig, ik).gt.0.) then
1383                              if (close_bottom(ig, ik).and.print_closure_warnings) then
1384                                    print *, "NOTICE: opening bottom of soil layer due to outward flux", ig, ik     
1385                              endif
1386                              close_bottom(ig, ik) = .false.
1387                        endif
1388                       
1389!                        if ( (flux(ig, ik).lt.flux(ig, ik - 1)).and.(flux(ig, ik).gt.0D0) ) then
1390!                              if (.not.close_top(ig, ik).and.print_closure_warnings) then
1391!                                    print *, "NOTICE: closing top of soil layer due to ice", ig, ik     
1392!                              endif
1393!                              close_top(ig, ik) = .true.
1394!                             
1395!                        elseif ( (flux(ig, ik).lt.flux(ig, ik - 1)).and.(flux(ig, ik - 1).lt.0D0) ) then
1396!                              if (.not.close_bottom(ig, ik).and.print_closure_warnings) then
1397!                                    print *, "NOTICE: closing bottom of soil layer due to ice", ig, ik
1398!                              endif
1399!                              close_bottom(ig, ik) = .true.
1400!                        endif
1401!                       
1402!                        if (close_top(ig, ik).and.close_bottom(ig, ik)) then
1403!                              close_top(ig, ik) = .false.
1404!                              close_bottom(ig, ik) = .false.
1405!                              if (print_closure_warnings) then
1406!                                    print *, "WARNING: Reopening soil layer after complete closure:", ig, ik
1407!                              endif
1408!                        elseif (close_top(ig, ik).or.close_bottom(ig, ik)) then
1409!                              if ((flux(ig, ik) - flux(ig, ik - 1)).gt.0D0) then
1410!                                    close_top(ig, ik) = .false.
1411!                                    close_bottom(ig, ik) = .false.
1412!                                    if (print_closure_warnings) then
1413!                                          print *, "WARNING: Reopening soil layer due to rising ice:", ig, ik
1414!                                    endif
1415!                              endif
1416!                             
1417!                              if (close_top(ig, ik).and.(flux(ig, ik - 1).lt.0D0)) then
1418!                                    close_top = .false.
1419!                                    if (print_closure_warnings) then
1420!                                          print *, "NOTICE: Reopening top of soil layer due to upward flux:", ig, ik
1421!                                    endif
1422!                              endif
1423!                             
1424!                              if (close_bottom(ig, ik).and.(flux(ig, ik).gt.0D0)) then
1425!                                    close_bottom = .false.
1426!                                    if (print_closure_warnings) then
1427!                                          print *, "NOTICE: Reopening bottom of soil layer due to downward flux:", ig, ik
1428!                                    endif
1429!                              endif
1430!                        endif
1431                  else  ! opening condition that ice content has fallen sufficiently
1432                        if (close_top(ig, ik).or.close_bottom(ig, ik).and.print_closure_warnings) then
1433                              print *, "NOTICE: Reopening soillayer due to decrease in ice", ig, ik                             
1434                        endif
1435                        close_top(ig, ik) = .false.
1436                        close_bottom(ig, ik) = .false.
1437                  endif
1438            enddo
1439      endif
1440enddo
1441
1442! Computation of total mass
1443! ===========================
1444
1445do ig = 1, ngrid
1446      ! initialize the masses to zero
1447      mass_h2o(ig) = 0.D0
1448      mass_ice(ig) = 0.D0
1449
1450      if(.not.watercaptag(ig)) then
1451            do ik = 1, nsoil
1452                  ! calculate the total water and ice mass
1453                  mass_h2o(ig) = mass_h2o(ig) +  h2otot(ig, ik) * interlayer_dz(ig, ik)
1454                  mass_ice(ig) = mass_ice(ig) +  ice(ig, ik) * interlayer_dz(ig, ik)
1455                 
1456                  ! calculate how close the water vapor content is to saturizing the adsorbed water
1457                  preduite(ig, ik) = znsoil(ig, ik) / nsat(ig, ik)
1458                 
1459                  ! write the results to the return variable
1460                  pqsoil(ig, ik, igcm_h2o_vap_soil) = znsoil(ig, ik)
1461                  pqsoil(ig, ik, igcm_h2o_ice_soil) = ice(ig, ik)
1462                  pqsoil(ig, ik, igcm_h2o_vap_ads) = adswater(ig, ik)
1463                 
1464
1465                  if (close_top(ig, ik)) then
1466                        close_out(ig, ik) = 1
1467                  elseif (close_bottom(ig, ik)) then
1468                        close_out(ig, ik) = -1
1469                  else
1470                        close_out(ig, ik) = 0 
1471                  endif
1472            enddo
1473      endif
1474     
1475      if (watercaptag(ig)) then
1476            do ik = 1, nsoil
1477                  saturation_water_ice(ig, ik) = -1
1478            enddo
1479      endif
1480     
1481      ! convert the logical :: flag exchange to a numeric output
1482      if (watercaptag(ig)) then
1483            exch(ig) = -2.
1484      elseif (exchange(ig)) then
1485            exch(ig) = 1.
1486      else
1487            exch(ig) = -1.
1488      endif
1489     
1490enddo
1491
1492! calculate the global total value
1493mass_h2o_tot = 0.D0
1494mass_ice_tot = 0.D0
1495do ig = 1, ngrid
1496      mass_h2o_tot = mass_h2o_tot + mass_h2o(ig) * cell_area(ig)
1497      mass_ice_tot = mass_ice_tot + mass_ice(ig) * cell_area(ig)
1498enddo
1499
1500! print and iterate the call number
1501! print * , 'Subsurface call n = ', n
1502n = n + 1
1503
1504! -----------------------------------------------------------------------
1505!      Output in diagfi and diagsoil
1506! -----------------------------------------------------------------------
1507
1508! reformat flux, because it has a unusual shape
1509do ig = 1, ngrid
1510      nsurf(ig) = rhosurf(ig) * pq(ig, 1, igcm_h2o_vap)
1511     
1512     
1513      do ik = 1, nsoil - 1
1514              var_flux_soil(ig, ik) = flux(ig, ik) 
1515      enddo
1516      var_flux_soil(ig, nsoil) = 0.
1517enddo
1518
1519if(writeoutput) then
1520!print *, "flux shape: ", shape(flux)
1521!print *, "adswater shape ", shape(adswater)
1522
1523!print *, "ngrid= ", ngrid
1524
1525!if (ngrid.eq.1) then  ! if in 1D
1526     
1527!      print *, "writing 1D data"
1528     
1529!      print *, real(adswater(:, :), 4)
1530!      print *, shape(adswater(:, :))
1531     
1532      zalpha(1, nsoil) = 0
1533      zbeta(1, nsoil) = 0
1534     
1535!       call write_output("zalpha","diffusion alpha coefficient", "unknown",real(zalpha(1, :)))
1536
1537!       call write_output("zbeta","diffusion beta coefficient", "unknown",real(zbeta(1, :)))
1538
1539!       call write_output("adswater","adswater", "kg / m^3",real(adswater(1, :)))
1540     
1541!       call write_output("znsoil","znsoil", "kg / m^3",real(znsoil(1, :)))
1542
1543!       call write_output("ice","ice", "kg / m^3",real(ice(1, :)))
1544     
1545!       call write_output("h2otot","total water content", "kg / m^3",real(h2otot(1, :)))
1546
1547!       call write_output("flux_soil","flux_soil", "kg / m^2",real(flux(1, :)))
1548     
1549!       call write_output("flux","flux soil", "kg / m^2",real(zdqsdifrego(:)))     
1550 
1551!       call write_output("flux_surf","surface flux", "kg / m^2",real(zdqsdifrego(1)))       
1552
1553!       call write_output("exchange","exchange", "boolean",real(exch(1)))               
1554
1555!else
1556     
1557!      print *, mass_h2o_tot
1558!      print *, real(mass_h2o_tot, 4)
1559     
1560!       call write_output("dztot1","change in dztot", "unknown",real(dztot1))
1561
1562        call write_output("flux_soillayer","flux of water between the soil layers", "kg m-2 s-1",var_flux_soil)                               
1563     
1564!       call write_output("A_coef","A coefficient of dztot1", "unknown",real(B))                               
1565     
1566!       call write_output("B_coef","B coefficient of dztot1", "unknown",real(B))       
1567
1568!       call write_output("H2O_init","initialized H2O", "kg/m^2",real(H2O))                                                         
1569
1570!      call write_output("H2O_depth_init","initialized H2O depth ", "m",real(H2O_depth))   
1571
1572       call write_output("ice_saturation_soil","Water ice saturation in the soil layers", "Percent",saturation_water_ice)
1573
1574!       call write_output("mass_h2o_tot","total mass of subsurface water over the planet", "kg",real(mass_h2o_tot))       
1575
1576!       call write_output("mass_ice_tot","total mass of subsurface ice over the planet", "kg",real(mass_ice_tot))         
1577
1578        call write_output("mass_h2o_soil","Mass of subsurface water column at each point", "kg m-2",(mass_h2o))                   
1579
1580        call write_output("mass_ice_soil","Mass of subsurface ice at each point", "kg m-2",(mass_ice))                   
1581           
1582        call write_output("znsoil","Water vapor soil concentration ", "kg.m - 3 of pore air",znsoil)                   
1583     
1584       call write_output("nsatsoil","subsurface water vapor saturation density", "kg/m^3",real(nsatsoil))                   
1585
1586       call write_output("nsurf","surface water vapor density", "kg/m^3",real(nsurf))                   
1587     
1588       call write_output("adswater","subsurface adsorbed water", "kg/m^3",real(adswater))                   
1589     
1590!       call write_output("subsurfice","subsurface ice", "kg/m^3",real(ice))                   
1591
1592        call write_output("flux_rego",'flux of water from the regolith','kg/m^2',zdqsdifrego)                   
1593     
1594!       call write_output("exchange",'exchange','no unit',real(exch))     
1595       
1596!       call write_output("close",'close top, bottom, or none (1, -1, 0)','no unit',real(close_out))                         
1597
1598!       call write_output("coeff_diffusion_soil",'interlayer diffusion coefficient','m^2/s',D_inter)                         
1599           
1600!endif
1601endif
1602RETURN
1603END
1604
1605
Note: See TracBrowser for help on using the repository browser.