1 | module 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 | USE comsoil_h_PEM, only: mlayer_PEM,layer_PEM ! Depth of the vertical grid |
---|
54 | USE soil_thermalproperties_mod, only: ice_thermal_properties |
---|
55 | implicit none |
---|
56 | ! inputs |
---|
57 | ! ----------- |
---|
58 | integer,intent(in) :: ngrid,nslope,nsoil_PEM ! Size of the physical grid, number of subslope, number of soil layer in the PEM |
---|
59 | logical,intent(in) :: watercaptag(ngrid) ! Boolean to check the presence of a perennial glacier |
---|
60 | real,intent(in) :: rhowatersurf_ave(ngrid,nslope) ! Water density at the surface, yearly averaged [kg/m^3] |
---|
61 | 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] |
---|
62 | real,intent(in) :: regolith_inertia(ngrid,nslope) ! Thermal inertia of the regolith layer [SI] |
---|
63 | ! Ouputs |
---|
64 | ! ----------- |
---|
65 | real,intent(out) :: ice_table_beg(ngrid,nslope) ! ice table depth [m] |
---|
66 | real,intent(out) :: ice_table_thickness(ngrid,nslope) ! ice table thickness [m] |
---|
67 | ! Local |
---|
68 | ! ----------- |
---|
69 | integer ig, islope,isoil,isoilend ! loop variables |
---|
70 | real :: diff_rho(nsoil_PEM) ! difference of water vapor density between the surface and at depth [kg/m^3] |
---|
71 | real :: ice_table_end ! depth of the end of the ice table [m] |
---|
72 | real :: previous_icetable_depth(ngrid,nslope) ! Ice table computed at previous ice depth [m] |
---|
73 | real :: stretch ! stretch factor to improve the convergence of the ice table |
---|
74 | real :: wice_inertia ! Water Ice thermal Inertia [USI] |
---|
75 | ! Code |
---|
76 | ! ----------- |
---|
77 | |
---|
78 | previous_icetable_depth(:,:) = ice_table_beg(:,:) |
---|
79 | do ig = 1,ngrid |
---|
80 | if(watercaptag(ig)) then |
---|
81 | ice_table_beg(ig,:) = 0. |
---|
82 | ice_table_thickness(ig,:) = layer_PEM(nsoil_PEM) ! Let's assume an infinite ice table (true when geothermal flux is set to 0.) |
---|
83 | else |
---|
84 | do islope = 1,nslope |
---|
85 | ice_table_beg(ig,islope) = -1. |
---|
86 | ice_table_thickness(ig,islope) = 0. |
---|
87 | do isoil = 1,nsoil_PEM |
---|
88 | diff_rho(isoil) = rhowatersurf_ave(ig,islope) - rhowatersoil_ave(ig,isoil,islope) |
---|
89 | enddo |
---|
90 | if(diff_rho(1) > 0) then ! ice is at the surface |
---|
91 | ice_table_beg(ig,islope) = 0. |
---|
92 | do isoilend = 2,nsoil_PEM-1 |
---|
93 | if((diff_rho(isoilend).gt.0).and.(diff_rho(isoilend+1).lt.0.)) then |
---|
94 | call findroot(diff_rho(isoilend),diff_rho(isoilend+1),mlayer_PEM(isoilend),mlayer_PEM(isoilend+1),ice_table_end) |
---|
95 | ice_table_thickness(ig,islope) = ice_table_end - ice_table_beg(ig,islope) |
---|
96 | exit |
---|
97 | endif |
---|
98 | enddo |
---|
99 | else |
---|
100 | 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 |
---|
101 | if((diff_rho(isoil).lt.0).and.(diff_rho(isoil+1).gt.0.)) then |
---|
102 | call findroot(diff_rho(isoil),diff_rho(isoil+1),mlayer_PEM(isoil),mlayer_PEM(isoil+1),ice_table_beg(ig,islope)) |
---|
103 | ! Now let's find the end of the ice table |
---|
104 | ice_table_thickness(ig,islope) = layer_PEM(nsoil_PEM) ! Let's assume an infinite ice table (true when geothermal flux is set to 0.) |
---|
105 | do isoilend = isoil+1,nsoil_PEM-1 |
---|
106 | if((diff_rho(isoilend).gt.0).and.(diff_rho(isoilend+1).lt.0.)) then |
---|
107 | call findroot(diff_rho(isoilend),diff_rho(isoilend+1),mlayer_PEM(isoilend),mlayer_PEM(isoilend+1),ice_table_end) |
---|
108 | ice_table_thickness(ig,islope) = ice_table_end - ice_table_beg(ig,islope) |
---|
109 | exit |
---|
110 | endif |
---|
111 | enddo |
---|
112 | endif !diff_rho(z) <0 & diff_rho(z+1) > 0 |
---|
113 | enddo |
---|
114 | endif ! diff_rho(1) > 0 |
---|
115 | enddo |
---|
116 | endif ! watercaptag |
---|
117 | enddo |
---|
118 | |
---|
119 | ! Small trick to speed up the convergence, Oded's idea. |
---|
120 | do islope = 1,nslope |
---|
121 | do ig = 1,ngrid |
---|
122 | if((ice_table_beg(ig,islope).gt.previous_icetable_depth(ig,islope)).and.(previous_icetable_depth(ig,islope).ge.0)) then |
---|
123 | call ice_thermal_properties(.false.,1.,regolith_inertia(ig,islope),wice_inertia) |
---|
124 | stretch = (regolith_inertia(ig,islope)/wice_inertia)**2 |
---|
125 | |
---|
126 | ice_table_thickness(ig,islope) = ice_table_thickness(ig,islope) + (ice_table_beg(ig,islope) - & |
---|
127 | previous_icetable_depth(ig,islope)+(ice_table_beg(ig,islope) - previous_icetable_depth(ig,islope))/stretch) |
---|
128 | ice_table_beg(ig,islope) = previous_icetable_depth(ig,islope)+(ice_table_beg(ig,islope) - previous_icetable_depth(ig,islope))/stretch |
---|
129 | endif |
---|
130 | enddo |
---|
131 | enddo |
---|
132 | |
---|
133 | RETURN |
---|
134 | END |
---|
135 | |
---|
136 | !!! -------------------------------------- |
---|
137 | |
---|
138 | |
---|
139 | SUBROUTINE compute_massh2o_exchange_ssi(ngrid,nslope,nsoil_PEM,former_ice_table_thickness,new_ice_table_thickness,ice_table_depth,tsoil,delta_m_h2o) |
---|
140 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
141 | !!! |
---|
142 | !!! Purpose: Compute the mass of H2O that has sublimated from the ice table / condensed |
---|
143 | !!! |
---|
144 | !!! Author: LL |
---|
145 | !!! |
---|
146 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
147 | use comsoil_h_PEM, only: mlayer_PEM |
---|
148 | use comslope_mod, only: subslope_dist,def_slope_mean |
---|
149 | use constants_marspem_mod,only:porosity |
---|
150 | #ifndef CPP_STD |
---|
151 | use comcstfi_h, only: pi |
---|
152 | #else |
---|
153 | use comcstfi_mod, only: pi |
---|
154 | #endif |
---|
155 | |
---|
156 | implicit none |
---|
157 | ! inputs |
---|
158 | ! ----------- |
---|
159 | integer,intent(in) :: ngrid,nslope,nsoil_PEM ! Size of the physical grid, number of subslope, number of soil layer in the PEM |
---|
160 | real,intent(in) :: former_ice_table_thickness(ngrid,nslope) ! ice table thickness at the former iteration [m] |
---|
161 | real,intent(in) :: new_ice_table_thickness(ngrid,nslope) ! ice table thickness at the current iteration [m] |
---|
162 | real,intent(out) :: ice_table_depth(ngrid,nslope) ! ice table depth [m] |
---|
163 | real,intent(in) :: tsoil(ngrid,nsoil_PEM,nslope) ! Soil temperature [K] |
---|
164 | ! outputs |
---|
165 | ! ----------- |
---|
166 | real,intent(out) :: delta_m_h2o(ngrid) ! Mass of H2O ice that has been condensed on the ice table / sublimates from the ice table [kg/m^2] |
---|
167 | |
---|
168 | ! local |
---|
169 | real :: rho(ngrid,nslope) ! density of water ice [kg/m^3] |
---|
170 | integer :: ig,islope,ilay,iref ! loop index |
---|
171 | real :: Tice(ngrid,nslope) ! ice temperature [k] |
---|
172 | ! Code |
---|
173 | ! ----------- |
---|
174 | rho(:,:) = 0. |
---|
175 | Tice(:,:) = 0. |
---|
176 | !1. First let's compute Tice using a linear interpolation between the mlayer level |
---|
177 | do ig = 1,ngrid |
---|
178 | do islope = 1,nslope |
---|
179 | if(ice_table_depth(ig,islope).gt.0.) then |
---|
180 | if (ice_table_depth(ig,islope).lt. 1e-10) then |
---|
181 | Tice(ig,islope) = tsoil(ig,1,islope) |
---|
182 | else |
---|
183 | iref=0 ! initialize iref |
---|
184 | do ilay=1,nsoil_PEM ! loop on layers to find the beginning of the ice table |
---|
185 | if (ice_table_depth(ig,islope).ge.mlayer_PEM(ilay)) then |
---|
186 | iref=ilay ! pure regolith layer up to here |
---|
187 | else |
---|
188 | ! correct iref was obtained in previous cycle |
---|
189 | exit |
---|
190 | endif |
---|
191 | enddo |
---|
192 | if(iref.eq.nsoil_PEM) then |
---|
193 | Tice(ig,islope) = tsoil(ig,nsoil_PEM,islope) |
---|
194 | else |
---|
195 | Tice(ig,islope) = (tsoil(ig,iref+1,islope) - tsoil(ig,iref,islope))/(mlayer_PEM(iref+1) - mlayer_PEM(iref))* & |
---|
196 | (ice_table_depth(ig,islope) - mlayer_PEM(iref)) + tsoil(ig,iref,islope) |
---|
197 | endif |
---|
198 | rho(ig,islope) = -3.5353e-4*Tice(ig,islope)**2+ 0.0351* Tice(ig,islope) + 933.5030 ! Rottgers, 2012 |
---|
199 | endif |
---|
200 | endif |
---|
201 | enddo |
---|
202 | enddo |
---|
203 | |
---|
204 | |
---|
205 | !2. Let's compute the amount of ice that has sublimated in each subslope |
---|
206 | do ig = 1,ngrid |
---|
207 | do islope = 1,nslope |
---|
208 | delta_m_h2o(ig) = delta_m_h2o(ig) + porosity*rho(ig,islope)*(new_ice_table_thickness(ig,islope) - former_ice_table_thickness(ig,islope)) & ! convention > 0. <=> it condenses |
---|
209 | *subslope_dist(ig,islope)/cos(def_slope_mean(islope)*pi/180.) |
---|
210 | enddo |
---|
211 | enddo |
---|
212 | |
---|
213 | return |
---|
214 | end subroutine |
---|
215 | |
---|
216 | !!! -------------------------------------- |
---|
217 | |
---|
218 | 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) |
---|
219 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
220 | !!! |
---|
221 | !!! Purpose: Compute layering between dry soil, pore filling ice, and ice sheet based on Schorgofer, Icarus (2010). Adapted from NS MSIM |
---|
222 | !!! |
---|
223 | !!! Author: LL |
---|
224 | !!! |
---|
225 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
226 | use comsoil_h_PEM,only: nsoilmx_PEM,mlayer_PEM |
---|
227 | use math_mod, only: deriv1,deriv1_onesided,colint,findroot,deriv2_simple |
---|
228 | implicit none |
---|
229 | ! inputs |
---|
230 | ! ------ |
---|
231 | |
---|
232 | real,intent(in) :: porefill(nsoilmx_PEM) ! Fraction of pore space filled with ice [Unitless] 0 <= f <= 1 for pore ice |
---|
233 | real,intent(in) :: psat_soil(nsoilmx_PEM) ! Soil water pressure at saturation, yearly averaged [Pa] |
---|
234 | real,intent(in) :: psat_surf ! surface water pressure at saturation, yearly averaged [Pa] |
---|
235 | real,intent(in) :: pwat_surf ! Water vapor pressure at the surface, not necesseraly at saturation, yearly averaged [Pa] |
---|
236 | real,intent(in) :: psat_bottom ! Boundary conditions for soil vapor pressure [Pa] |
---|
237 | real,intent(in) :: B ! constant (Eq. 8 from Schorgofer, Icarus (2010).) |
---|
238 | integer, intent(in) :: index_IS ! index of the soil layer where the ice sheet begins [1] |
---|
239 | real, intent(inout) :: depth_filling ! depth where pore filling begins [m] |
---|
240 | |
---|
241 | ! outputs |
---|
242 | ! ------- |
---|
243 | integer,intent(out) :: index_filling ! index where the pore filling begins [1] |
---|
244 | integer, intent(out) :: index_geothermal ! index where the ice table stops [1] |
---|
245 | real, intent(out) :: depth_geothermal ! depth where the ice table stops [m] |
---|
246 | real, intent(out) :: dz_etadz_rho(nsoilmx_PEM) ! \partial z(eta \partial z rho), eta is the constriction, used later for pore filling increase |
---|
247 | |
---|
248 | ! local |
---|
249 | |
---|
250 | real :: eta(nsoilmx_PEM) ! constriction |
---|
251 | integer :: ilay ! index for loop |
---|
252 | real :: old_depth_filling ! depth_filling saved |
---|
253 | real :: dz_psat(nsoilmx_PEM) ! first derivative of the vapor pressure at saturationn |
---|
254 | integer :: index_tmp ! for loop |
---|
255 | real :: Jdry ! flux trought the dry layer |
---|
256 | real :: Jsat ! flux trought the ice layer |
---|
257 | real :: Jdry_prevlay,Jsat_prevlay ! same but for the previous ice layer |
---|
258 | integer :: index_firstice ! first index where ice appears (i.e., f > 0) |
---|
259 | real :: dz_eta(nsoilmx_PEM) ! \partial z \eta |
---|
260 | real :: dzz_psat(nsoilmx_PEM) ! \partial \partial psat |
---|
261 | real :: massfillabove,massfillafter ! h2O mass above and after index_geothermal |
---|
262 | |
---|
263 | ! constant |
---|
264 | real :: pvap2rho = 18.e-3/8.314 |
---|
265 | ! |
---|
266 | |
---|
267 | ! 0. Compute constriction over the layer |
---|
268 | ! Within the ice sheet, constriction is set to 0. Elsewhere, constriction = (1-porefilling)**2 |
---|
269 | if (index_IS.lt.0) then |
---|
270 | index_tmp = nsoilmx_PEM |
---|
271 | do ilay = 1,nsoilmx_PEM |
---|
272 | call constriction(porefill(ilay),eta(ilay)) |
---|
273 | enddo |
---|
274 | else |
---|
275 | index_tmp = index_IS |
---|
276 | do ilay = 1,index_IS-1 |
---|
277 | call constriction(porefill(ilay),eta(ilay)) |
---|
278 | enddo |
---|
279 | do ilay = index_IS,nsoilmx_PEM |
---|
280 | eta(ilay) = 0. |
---|
281 | enddo |
---|
282 | endif |
---|
283 | |
---|
284 | ! 1. Depth at which pore filling occurs. We solve Eq. 9 from Schorgofer, Icarus (2010) |
---|
285 | |
---|
286 | old_depth_filling = depth_filling |
---|
287 | |
---|
288 | call deriv1(mlayer_PEM,nsoilmx_PEM,psat_soil,psat_surf,psat_bottom,dz_psat) |
---|
289 | |
---|
290 | do ilay = 1,index_tmp |
---|
291 | Jdry = (psat_soil(ilay) - pwat_surf)/mlayer_PEM(ilay) ! left member of Eq. 9 from Schorgofer, Icarus (2010) |
---|
292 | Jsat = eta(ilay)*dz_psat(ilay) !right member of Eq. 9 from Schorgofer, Icarus (2010) |
---|
293 | if((Jdry - Jsat).le.0) then |
---|
294 | index_filling = ilay |
---|
295 | exit |
---|
296 | endif |
---|
297 | enddo |
---|
298 | |
---|
299 | if(index_filling.eq.1) depth_filling = mlayer_PEM(1) |
---|
300 | |
---|
301 | if(index_filling.gt.1) then |
---|
302 | Jdry_prevlay = (psat_soil(index_filling-1) - pwat_surf)/mlayer_PEM(index_filling-1) |
---|
303 | Jsat_prevlay = eta(index_filling-1)*dz_psat(index_filling-1) |
---|
304 | call findroot(Jdry-Jsat,Jdry_prevlay-Jsat_prevlay,mlayer_PEM(index_filling),mlayer_PEM(index_filling-1),depth_filling) |
---|
305 | endif |
---|
306 | |
---|
307 | |
---|
308 | ! 2. Compute d_z (eta* d_z(rho)) (last term in Eq. 13 of Schorgofer, Icarus (2010)) |
---|
309 | |
---|
310 | |
---|
311 | ! 2.0 preliminary: depth to shallowest ice (discontinuity at interface) |
---|
312 | |
---|
313 | index_firstice = -1 |
---|
314 | do ilay = 1,nsoilmx_PEM |
---|
315 | if (porefill(ilay).le.0.) then |
---|
316 | index_firstice = ilay ! first point with ice |
---|
317 | exit |
---|
318 | endif |
---|
319 | enddo |
---|
320 | |
---|
321 | ! 2.1: now we can computeCompute d_z (eta* d_z(rho)) |
---|
322 | |
---|
323 | call deriv1(mlayer_PEM,nsoilmx_PEM,eta,1.,eta(nsoilmx_PEM-1),dz_eta) |
---|
324 | |
---|
325 | if ((index_firstice.gt.0).and.(index_firstice.lt.nsoilmx_PEM-2)) then |
---|
326 | call deriv1_onesided(index_firstice,mlayer_PEM,nsoilmx_PEM,eta,dz_eta(index_firstice)) |
---|
327 | endif |
---|
328 | |
---|
329 | call deriv2_simple(mlayer_PEM,nsoilmx_PEM,psat_soil,psat_surf,psat_bottom,dzz_psat) |
---|
330 | dz_etadz_rho(:) = pvap2rho*(dz_eta(:)*dz_psat(:) + eta(:)*dzz_psat(:)) |
---|
331 | |
---|
332 | ! 3. Ice table boundary due to geothermal heating |
---|
333 | |
---|
334 | if(index_IS.gt.0) index_geothermal = -1 |
---|
335 | if(index_geothermal.lt.0) depth_geothermal = -1. |
---|
336 | if((index_geothermal.gt.0).and.(index_IS.lt.0)) then ! Eq. 21 from Schorfoger, Icarus (2010) |
---|
337 | index_geothermal = -1 |
---|
338 | do ilay=2,nsoilmx_PEM |
---|
339 | if (dz_psat(ilay).gt.0.) then ! first point with reversed flux |
---|
340 | index_geothermal=ilay |
---|
341 | call findroot(dz_psat(ilay-1),dz_psat(ilay),mlayer_PEM(ilay-1),mlayer_PEM(ilay),depth_geothermal) |
---|
342 | exit |
---|
343 | endif |
---|
344 | enddo |
---|
345 | else |
---|
346 | index_geothermal = -1 |
---|
347 | endif |
---|
348 | if ((index_geothermal.gt.0).and.(index_IS.lt.0)) then ! Eq. 24 from Schorgofer, Icarus (2010) |
---|
349 | call colint(porefill(:)/eta(:),mlayer_PEM,nsoilmx_PEM,index_geothermal-1,nsoilmx_PEM,massfillabove) |
---|
350 | index_tmp = -1 |
---|
351 | do ilay=index_geothermal,nsoilmx_PEM |
---|
352 | if (minval(eta(ilay:nsoilmx_PEM)).le.0.) cycle ! eta=0 means completely full |
---|
353 | call colint(porefill(:)/eta(:),mlayer_PEM,nsoilmx_PEM,ilay,nsoilmx_PEM,massfillafter) |
---|
354 | if (massfillafter<dz_psat(ilay)*pvap2rho*B) then ! usually executes on i=typeG |
---|
355 | if (ilay>index_geothermal) then |
---|
356 | ! write(34,*) '# adjustment to geotherm depth by',ilay-index_geothermal |
---|
357 | call findroot(dz_psat(ilay-1)*pvap2rho*B-massfillabove, & |
---|
358 | dz_psat(ilay)*pvap2rho*B-massfillafter,mlayer_PEM(ilay-1),mlayer_PEM(ilay),depth_geothermal) |
---|
359 | ! if (depth_geothermal.gt.mlayer_PEM(ilay) .or. depth_geothermal.lt.<mlayer_PEM(ilay-1)) write(34,*) |
---|
360 | ! '# WARNING: zdepthG interpolation failed',ilay,mlayer_PEM(ilay-1),depth_geothermal,mlayer_PEM(ilay) |
---|
361 | index_tmp=ilay |
---|
362 | endif |
---|
363 | ! otherwise leave depth_geothermal unchanged |
---|
364 | exit |
---|
365 | endif |
---|
366 | massfillabove = massfillafter |
---|
367 | enddo |
---|
368 | if (index_tmp.gt.0) index_geothermal = index_tmp |
---|
369 | end if |
---|
370 | return |
---|
371 | end subroutine |
---|
372 | |
---|
373 | !!! -------------------------------------- |
---|
374 | |
---|
375 | SUBROUTINE constriction(porefill,eta) |
---|
376 | |
---|
377 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
378 | !!! |
---|
379 | !!! Purpose: Compute the constriction of vapor flux by pore ice |
---|
380 | !!! |
---|
381 | !!! Author: LL |
---|
382 | !!! |
---|
383 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
384 | implicit none |
---|
385 | real,intent(in) :: porefill ! pore filling fraction |
---|
386 | real,intent(out) :: eta ! constriction |
---|
387 | |
---|
388 | !!! |
---|
389 | |
---|
390 | if (porefill.le.0.) eta = 1. |
---|
391 | if ((porefill.gt.0.) .and.(porefill.lt.1.)) then |
---|
392 | eta = (1-porefill)**2 ! Hudson et al., JGR, 2009 |
---|
393 | endif |
---|
394 | if (porefill.le.1.) eta = 0. |
---|
395 | return |
---|
396 | end subroutine |
---|
397 | |
---|
398 | end module |
---|