source: trunk/LMDZ.COMMON/libf/evolution/ice_table_mod.F90 @ 2982

Last change on this file since 2982 was 2980, checked in by romain.vande, 3 years ago

Mars PEM :

Adapt PEM to 1d runs.
Cleaning of names and unused variables.
Correct minor errors.
Adapt and correct reshape_xios_output utilitary for 1d diagfi output.

RV

File size: 14.4 KB
Line 
1module ice_table_mod
2
3  implicit none
4
5!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
6!!!
7!!! Purpose: Ice table (pore-filling) variables and methods to compute it (dynamic and static)
8!!! Author:  LL, 02/2023
9!!!
10!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
11
12  LOGICAL,SAVE :: icetable_equilibrium                    ! Boolean to say if the PEM needs to recompute the icetable depth when at  equilibrium
13  LOGICAL,SAVE :: icetable_dynamic                        ! Boolean to say if the PEM needs to recompute the icetable depth (dynamic method)
14  real,save,allocatable  :: porefillingice_depth(:,:)    ! ngrid x nslope: Depth of the ice table [m]
15  real,save,allocatable  :: porefillingice_thickness(:,:)    ! ngrid x nslope: Thickness of the ice table [m]
16
17  contains
18
19
20  subroutine ini_ice_table_porefilling(ngrid,nslope)
21 
22    implicit none
23    integer,intent(in) :: ngrid ! number of atmospheric columns
24    integer,intent(in) :: nslope ! number of slope within a mesh
25
26    allocate(porefillingice_depth(ngrid,nslope))
27    allocate(porefillingice_thickness(ngrid,nslope))
28
29  end subroutine ini_ice_table_porefilling
30
31  subroutine end_ice_table_porefilling
32
33    implicit none
34    if (allocated(porefillingice_depth)) deallocate(porefillingice_depth)
35    if (allocated(porefillingice_thickness)) deallocate(porefillingice_thickness)
36
37  end subroutine end_ice_table_porefilling
38
39!!! --------------------------------------
40
41   SUBROUTINE computeice_table_equilibrium(ngrid,nslope,nsoil_PEM,watercaptag,rhowatersurf_ave,rhowatersoil_ave,regolith_inertia,ice_table_beg,ice_table_thickness)
42!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
43!!!
44!!! Purpose: Compute the ice table depth (pore-filling) knowing the yearly average water
45!!! density at the surface and at depth.
46!!! Computations are made following the methods in Schorgofer et al., 2005
47!!! This subroutine only gives the ice table at equilibrium and does not consider exchange with the atmosphere
48!!!
49!!! Author: LL
50!!!
51!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
52    use math_mod,only: findroot
53#ifndef CPP_STD
54    USE comsoil_h_PEM, only: mlayer_PEM,layer_PEM                  ! Depth of the vertical grid
55    USE soil_thermalproperties_mod, only: ice_thermal_properties
56    implicit none
57!   inputs
58! -----------
59    integer,intent(in) :: ngrid,nslope,nsoil_PEM                   ! Size of the physical grid, number of subslope, number of soil layer in the PEM
60    logical,intent(in) :: watercaptag(ngrid)                       ! Boolean to check the presence of a perennial glacier
61    real,intent(in) :: rhowatersurf_ave(ngrid,nslope)              ! Water density at the surface, yearly averaged [kg/m^3]
62    real,intent(in) :: rhowatersoil_ave(ngrid,nsoil_PEM,nslope)    ! Water density at depth, computed from clapeyron law's (Murchy and Koop 2005), yearly averaged  [kg/m^3]
63    real,intent(in) :: regolith_inertia(ngrid,nslope)              ! Thermal inertia of the regolith layer [SI]
64!   Ouputs
65! -----------
66    real,intent(out) :: ice_table_beg(ngrid,nslope)               ! ice table depth [m]
67    real,intent(out) :: ice_table_thickness(ngrid,nslope)         ! ice table thickness [m]
68!   Local
69! -----------
70    real :: z1,z2                                                  ! intermediate variables used when doing a linear interpolation between two depths to find the root
71    integer ig, islope,isoil,isoilend                              ! loop variables
72    real :: diff_rho(nsoil_PEM)                                    ! difference of water vapor density between the surface and at depth [kg/m^3]
73    real  :: ice_table_end                                         ! depth of the end of the ice table  [m]
74    real :: previous_icetable_depth(ngrid,nslope)                  ! Ice table computed at previous ice depth [m]
75    real :: stretch                                                ! stretch factor to improve the convergence of the ice table
76    real :: wice_inertia                                           ! Water Ice thermal Inertia [USI]   
77!   Code
78! -----------
79
80     previous_icetable_depth(:,:) = ice_table_beg(:,:)
81     do ig = 1,ngrid
82      if(watercaptag(ig)) then
83         ice_table_beg(ig,:) = 0.
84         ice_table_thickness(ig,:) = layer_PEM(nsoil_PEM) ! Let's assume an infinite ice table (true when geothermal flux is set to 0.)
85      else
86        do islope = 1,nslope
87           ice_table_beg(ig,islope) = -1.
88           ice_table_thickness(ig,islope) = 0.
89         do isoil = 1,nsoil_PEM
90           diff_rho(isoil) = rhowatersurf_ave(ig,islope) - rhowatersoil_ave(ig,isoil,islope)
91         enddo
92         if(diff_rho(1) > 0) then                                   ! ice is at the surface
93           ice_table_beg(ig,islope) = 0.
94               do isoilend = 2,nsoil_PEM-1
95                 if((diff_rho(isoilend).gt.0).and.(diff_rho(isoilend+1).lt.0.)) then
96                  call findroot(diff_rho(isoilend),diff_rho(isoilend+1),mlayer_PEM(isoilend),mlayer_PEM(isoilend+1),ice_table_end)
97                  ice_table_thickness(ig,islope) = ice_table_end - ice_table_beg(ig,islope)
98                  exit
99                 endif
100               enddo
101         else
102           do isoil = 1,nsoil_PEM -1                                ! general case, we find the ice table depth by doing a linear approximation between the two depth, and then solve the first degree equation to find the root
103             if((diff_rho(isoil).lt.0).and.(diff_rho(isoil+1).gt.0.)) then
104               call findroot(diff_rho(isoil),diff_rho(isoil+1),mlayer_PEM(isoil),mlayer_PEM(isoil+1),ice_table_beg(ig,islope))
105               ! Now let's find the end of the ice table
106               ice_table_thickness(ig,islope) = layer_PEM(nsoil_PEM) ! Let's assume an infinite ice table (true when geothermal flux is set to 0.)
107               do isoilend = isoil+1,nsoil_PEM-1
108                 if((diff_rho(isoilend).gt.0).and.(diff_rho(isoilend+1).lt.0.)) then
109                  call findroot(diff_rho(isoilend),diff_rho(isoilend+1),mlayer_PEM(isoilend),mlayer_PEM(isoilend+1),ice_table_end)
110                  ice_table_thickness(ig,islope) = ice_table_end - ice_table_beg(ig,islope)
111                  exit
112                 endif
113               enddo
114             endif !diff_rho(z) <0 & diff_rho(z+1) > 0
115           enddo
116          endif ! diff_rho(1) > 0
117        enddo
118       endif ! watercaptag
119      enddo
120
121! Small trick to speed up the convergence, Oded's idea.
122      do islope = 1,nslope
123        do ig = 1,ngrid
124         if((ice_table_beg(ig,islope).gt.previous_icetable_depth(ig,islope)).and.(previous_icetable_depth(ig,islope).ge.0)) then
125            call ice_thermal_properties(.false.,1.,regolith_inertia(ig,islope),wice_inertia)
126            stretch = (regolith_inertia(ig,islope)/wice_inertia)**2
127
128            ice_table_thickness(ig,islope) = ice_table_thickness(ig,islope) + (ice_table_beg(ig,islope) - &
129            previous_icetable_depth(ig,islope)+(ice_table_beg(ig,islope) - previous_icetable_depth(ig,islope))/stretch)
130            ice_table_beg(ig,islope) = previous_icetable_depth(ig,islope)+(ice_table_beg(ig,islope) - previous_icetable_depth(ig,islope))/stretch
131         endif
132        enddo
133      enddo
134
135      RETURN
136#endif
137      END
138
139!!! --------------------------------------
140
141   SUBROUTINE find_layering_icetable(porefill,psat_soil,psat_surf,pwat_surf,psat_bottom,B,index_IS,depth_filling,index_filling,index_geothermal,depth_geothermal,dz_etadz_rho)
142!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
143!!!
144!!! Purpose: Compute layering between dry soil, pore filling ice, and ice sheet based on Schorgofer, Icarus (2010). Adapted from NS MSIM
145!!!
146!!! Author: LL
147!!!
148!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
149use comsoil_h_PEM,only: nsoilmx_PEM,mlayer_PEM
150use math_mod, only: deriv1,deriv1_onesided,colint,findroot,deriv2_simple
151implicit none
152! inputs
153! ------
154
155real,intent(in) :: porefill(nsoilmx_PEM) ! Fraction of pore space filled with ice [Unitless] 0 <= f <= 1 for pore ice
156real,intent(in) :: psat_soil(nsoilmx_PEM) ! Soil water pressure at saturation, yearly averaged [Pa]
157real,intent(in) :: psat_surf ! surface water pressure at saturation, yearly averaged [Pa]
158real,intent(in) :: pwat_surf ! Water vapor pressure  at the surface, not necesseraly at saturation, yearly averaged [Pa]
159real,intent(in) :: psat_bottom ! Boundary conditions for soil vapor pressure [Pa]
160real,intent(in) :: B ! constant (Eq. 8 from  Schorgofer, Icarus (2010).)
161integer, intent(in) :: index_IS ! index of the soil layer where the ice sheet begins [1]
162real, intent(inout) :: depth_filling ! depth where pore filling begins [m]
163
164! outputs
165! -------
166integer,intent(out) :: index_filling ! index where the pore filling begins [1]
167integer, intent(out) :: index_geothermal ! index where the ice table stops [1]
168real, intent(out) :: depth_geothermal ! depth where the ice table stops [m]
169real, intent(out) :: dz_etadz_rho(nsoilmx_PEM) ! \partial z(eta \partial z rho), eta is the constriction, used later for pore filling increase
170
171! local
172
173real :: eta(nsoilmx_PEM)  ! constriction
174integer :: ilay ! index for loop
175real :: old_depth_filling ! depth_filling saved
176real :: dz_psat(nsoilmx_PEM) ! first derivative of the vapor pressure at saturationn
177integer :: index_tmp ! for loop
178real :: Jdry ! flux trought the dry layer
179real :: Jsat ! flux trought the ice layer
180real :: Jdry_prevlay,Jsat_prevlay ! same but for the previous ice layer
181integer :: index_firstice ! first index where ice appears (i.e., f > 0)
182real :: dz_eta(nsoilmx_PEM) ! \partial z \eta
183real :: dz_eta_low ! same but evaluated at the interface for ice
184real :: dzz_psat(nsoilmx_PEM) ! \partial \partial psat
185real :: massfillabove,massfillafter ! h2O mass above and after index_geothermal
186
187! constant
188real :: pvap2rho = 18.e-3/8.314
189!
190 
191! 0. Compute constriction over the layer
192! Within the ice sheet, constriction is set to 0. Elsewhere, constriction =  (1-porefilling)**2
193if (index_IS.lt.0) then
194   index_tmp = nsoilmx_PEM
195   do ilay = 1,nsoilmx_PEM
196       call constriction(porefill(ilay),eta(ilay))
197   enddo
198else
199   index_tmp = index_IS
200   do ilay = 1,index_IS-1
201       call constriction(porefill(ilay),eta(ilay))
202   enddo
203   do ilay = index_IS,nsoilmx_PEM
204       eta(ilay) = 0.
205   enddo
206endif
207
208! 1. Depth at which pore filling occurs. We solve Eq. 9 from  Schorgofer, Icarus  (2010)   
209
210old_depth_filling = depth_filling
211
212call deriv1(mlayer_PEM,nsoilmx_PEM,psat_soil,psat_surf,psat_bottom,dz_psat) 
213
214do ilay = 1,index_tmp
215  Jdry = (psat_soil(ilay) - pwat_surf)/mlayer_PEM(ilay) ! left member of Eq. 9 from Schorgofer, Icarus  (2010)
216  Jsat =  eta(ilay)*dz_psat(ilay) !right member of Eq. 9 from Schorgofer, Icarus (2010)
217  if((Jdry - Jsat).le.0) then
218     index_filling = ilay
219     exit
220   endif
221enddo
222
223if(index_filling.eq.1)  depth_filling = mlayer_PEM(1)
224
225if(index_filling.gt.1) then
226    Jdry_prevlay = (psat_soil(index_filling-1) - pwat_surf)/mlayer_PEM(index_filling-1)
227    Jsat_prevlay = eta(index_filling-1)*dz_psat(index_filling-1)
228   call findroot(Jdry-Jsat,Jdry_prevlay-Jsat_prevlay,mlayer_PEM(index_filling),mlayer_PEM(index_filling-1),depth_filling)
229endif
230
231
232! 2. Compute d_z (eta* d_z(rho)) (last term in Eq. 13 of Schorgofer, Icarus (2010))
233
234
235! 2.0 preliminary: depth to shallowest  ice (discontinuity at interface)
236
237index_firstice = -1
238do ilay = 1,nsoilmx_PEM
239   if (porefill(ilay).le.0.) then
240        index_firstice = ilay  ! first point with ice
241        exit
242    endif
243enddo
244
245! 2.1: now we can computeCompute d_z (eta* d_z(rho))
246
247call deriv1(mlayer_PEM,nsoilmx_PEM,eta,1.,eta(nsoilmx_PEM-1),dz_eta)
248
249if ((index_firstice.gt.0).and.(index_firstice.lt.nsoilmx_PEM-2)) then
250     call deriv1_onesided(index_firstice,mlayer_PEM,nsoilmx_PEM,eta,dz_eta(index_firstice))
251endif
252
253call deriv2_simple(mlayer_PEM,nsoilmx_PEM,psat_soil,psat_surf,psat_bottom,dzz_psat)
254dz_etadz_rho(:) = pvap2rho*(dz_eta(:)*dz_psat(:) + eta(:)*dzz_psat(:))
255
256! 3. Ice table boundary due to geothermal heating
257
258if(index_IS.gt.0) index_geothermal = -1
259if(index_geothermal.lt.0) depth_geothermal = -1.
260if((index_geothermal.gt.0).and.(index_IS.lt.0)) then ! Eq. 21 from Schorfoger, Icarus (2010)
261   index_geothermal = -1
262   do ilay=2,nsoilmx_PEM
263        if (dz_psat(ilay).gt.0.) then  ! first point with reversed flux
264           index_geothermal=ilay
265           call findroot(dz_psat(ilay-1),dz_psat(ilay),mlayer_PEM(ilay-1),mlayer_PEM(ilay),depth_geothermal)
266           exit
267        endif
268     enddo
269  else
270     index_geothermal = -1
271  endif
272 if ((index_geothermal.gt.0).and.(index_IS.lt.0)) then !  Eq. 24 from Schorgofer, Icarus (2010)
273     call  colint(porefill(:)/eta(:),mlayer_PEM,nsoilmx_PEM,index_geothermal-1,nsoilmx_PEM,massfillabove)
274     index_tmp = -1
275     do ilay=index_geothermal,nsoilmx_PEM
276        if (minval(eta(ilay:nsoilmx_PEM)).le.0.) cycle  ! eta=0 means completely full
277         call colint(porefill(:)/eta(:),mlayer_PEM,nsoilmx_PEM,ilay,nsoilmx_PEM,massfillafter)
278        if (massfillafter<dz_psat(ilay)*pvap2rho*B) then  ! usually executes on i=typeG
279           if (ilay>index_geothermal) then
280!              write(34,*) '# adjustment to geotherm depth by',ilay-index_geothermal
281              call findroot(dz_psat(ilay-1)*pvap2rho*B-massfillabove, & 
282                   dz_psat(ilay)*pvap2rho*B-massfillafter,mlayer_PEM(ilay-1),mlayer_PEM(ilay),depth_geothermal)
283!               if (depth_geothermal.gt.mlayer_PEM(ilay) .or. depth_geothermal.lt.<mlayer_PEM(ilay-1)) write(34,*)
284!                     '# WARNING: zdepthG interpolation failed',ilay,mlayer_PEM(ilay-1),depth_geothermal,mlayer_PEM(ilay)
285              index_tmp=ilay
286           endif
287           ! otherwise leave depth_geothermal unchanged
288           exit
289        endif
290        massfillabove = massfillafter
291     enddo
292     if (index_tmp.gt.0) index_geothermal = index_tmp
293  end if
294  return
295end subroutine
296
297!!! --------------------------------------
298
299SUBROUTINE constriction(porefill,eta)
300
301        !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
302        !!!
303        !!! Purpose: Compute the  constriction of vapor flux by pore ice
304        !!!
305        !!! Author: LL
306        !!!
307        !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
308        implicit none
309        real,intent(in) :: porefill ! pore filling fraction
310        real,intent(out) :: eta ! constriction
311
312        !!!
313
314        if (porefill.le.0.) eta = 1.
315        if ((porefill.gt.0.) .and.(porefill.lt.1.)) then
316        eta = (1-porefill)**2  ! Hudson et al., JGR, 2009
317        endif
318        if (porefill.le.1.) eta = 0.
319        return
320        end subroutine
321
322end module
Note: See TracBrowser for help on using the repository browser.