source: trunk/LMDZ.COMMON/libf/evolution/subsurface_ice.F90 @ 4076

Last change on this file since 4076 was 4074, checked in by jbclement, 11 days ago

PEM:

  • Correct management of H2O ice tendency in 1D when there is not enough ice anymore.
  • Clean initialization of allocatable module arrays (especially needed when no slope)
  • One more renaming for consistency + few small updates thoughout the code.

JBC

File size: 30.9 KB
Line 
1MODULE subsurface_ice
2!-----------------------------------------------------------------------
3! NAME
4!     subsurface_ice
5!
6! DESCRIPTION
7!     Retreat and growth of subsurface ice on Mars with constant orbital
8!     elements.
9!
10! AUTHORS & DATE
11!     N. Schorghofer (MSIM dynamical program)
12!     E. Vos, 2025
13!     JB Clement, 2025
14!
15! NOTES
16!    Based on Norbert Schorgofer's code. Updated for PEM.
17!-----------------------------------------------------------------------
18
19! DEPENDENCIES
20! ------------
21use numerics, only: dp, di
22
23! DECLARATION
24! -----------
25implicit none
26
27! PARAMETERS
28! ----------
29real(8), parameter :: dt = 0.02  ! in units of Mars solar days
30!real(8), parameter :: Fgeotherm = 0.
31real(8), parameter :: Fgeotherm = 0 !0.028  ! [W/m^2]
32real(8), parameter :: Lco2frost=6.0e5, co2albedo=0.60, co2emiss=1.
33real(8), parameter :: emiss0 = 1.  ! emissivity of dry surface
34integer, parameter :: EQUILTIME =15 ! [Mars years]
35
36integer, parameter :: NMAX = 1000
37
38contains
39!++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
40
41
42!=======================================================================
43!============================ dyn_ss_ice_m =============================
44!=======================================================================
45
46
47SUBROUTINE dyn_ss_ice_m(ssi_depth_in,T1,Tb,nz,thIn,p0,pfrost,porefill_in,porefill,ssi_depth)
48
49!***********************************************************************
50! Retreat and growth of subsurface ice on Mars
51! orbital elements remain constant
52!***********************************************************************
53  use orbit, only: sol_len_s
54  use maths, only: pi
55  ! DECLARATION
56! -----------
57implicit none
58  integer, parameter :: NP=1   ! # of sites
59  integer nz, i, k, iloop
60  real(8)  zmax, delta, z(NMAX), icetime, porosity, icefrac
61  real(8), dimension(NP) ::  albedo, thIn, rhoc
62  real(8), dimension(NP) :: pfrost, p0
63  real(8) newti, stretch, newrhoc, ecc, omega, eps, timestep
64  real(8)  ssi_depth_in, ssi_depth, T1
65  real(8), dimension(NP) :: zdepthF, zdepthE, zdepthT, zdepthG
66  real(8), dimension(nz,NP) :: porefill, porefill_in
67  real(8), dimension(nz) :: Tb
68  real(8), dimension(NP) :: Tmean1, Tmean3, avrho1
69  real(8) tmax, tlast, avrho1prescribed(NP), l1
70  real(8), external :: smartzfac
71
72  !if (iargc() /= 1) then
73  !   stop 'USAGE: icages ext'
74  !endif
75  !call getarg( 1, ext )
76
77  if (NP>100) stop 'subroutine icelayer_mars cannot handle this many sites'
78
79  ! parameters that never ever change
80  porosity = 0.4d0  ! porosity of till
81  !rhoc(:) = 1500.*800.  ! will be overwritten
82  icefrac = 0.98
83  tmax = 1
84  tlast = 0.
85  avrho1prescribed(:) = pfrost/T1  ! <0 means absent
86  albedo=0.23
87  !avrho1prescribed(:) = 0.16/200.  ! units are Pa/K
88
89  !open(unit=21,file='lats.'//ext,action='read',status='old',iostat=ierr)
90  !if (ierr /= 0) then
91  !   print *,'File lats.'//ext,'not found'
92  !   stop
93  !endif
94  do k=1,NP
95     !read(21,*) latitude(k),albedo(k),thIn(k),htopo(k)
96     ! empirical relation from Mellon & Jakosky
97     rhoc(k) = 800.*(150.+100.*sqrt(34.2+0.714*thIn(k)))
98  enddo
99  !close(21)
100
101  ! set eternal grid
102  zmax = 25.
103  !zfac = smartzfac(nz,zmax,6,0.032d0)
104  !call setgrid(nz,z,zmax,zfac)
105  l1=2.e-4
106  do iloop=0,nz - 1
107    z(iloop + 1) = l1*(1+iloop**2.9*(1-exp(-real(iloop)/20.)))
108  enddo
109
110
111  !open(unit=30,file='z.'//ext,action='write',status='unknown')
112  !write(30,'(999(f8.5,1x))') z(1:nz)
113  !close(30)
114
115  !ecc = ecc_in;  eps = obl_in*d2r;  omega = Lp_in*d2r   ! today
116  ! total atmospheric pressure
117  !p0(:) = 600.
118  ! presently 520 Pa at zero elevation (Smith & Zuber, 1998)
119 ! do k=1,NP
120 !    p0(k)=520*exp(-htopo(k)/10800.)
121 ! enddo
122  timestep = 1 ! must be integer fraction of 1 ka
123  icetime = -tmax-timestep  ! earth years
124
125  ! initializations
126  !Tb = -9999.
127  zdepthF(:) = -9999.
128
129  !zdepthT(1:NP) = -9999.   ! reset again below
130 ! zdepthT(1:NP) = 0.
131
132 ! print *,'RUNNING MARS_FAST'
133 ! print *,'Global model parameters:'
134 ! print *,'nz=',nz,' zfac=',zfac,'zmax=',zmax
135 ! print *,'porosity=',porosity
136 ! print *,'starting at time',icetime,'years'
137 ! print *,'time step=',timestep,'years'
138 ! print *,'eps=',eps/d2r,'ecc=',ecc,'omega=',omega/d2r
139 ! print *,'number of sites=',NP
140 ! print *,'Site specific parameters:'
141  do k=1,NP
142     if (NP>1) print *,'  Site ',k
143 !    print *,'  latitude (deg)',latitude(k),' rho*c (J/m^3/K)',rhoc(k),' thIn=',thIn(k)
144 !    print *,'  total pressure=',p0(k),'partial pressure=',pfrost(k)
145     delta = thIn(k)/rhoc(k)*sqrt(sol_len_s/pi)
146 !    print *,'  skin depths (m)',delta,delta*sqrt(solsperyear)
147     call soilthprop(porosity,1.d0,rhoc(k),thIn(k),1,newrhoc,newti,icefrac)
148     stretch = (newti/thIn(k))*(rhoc(k)/newrhoc)
149     do i=1,nz
150        if (z(i)<delta) cycle
151 !       print *,'  ',i-1,' grid points within diurnal skin depth'
152        exit
153     enddo
154 !    print *,'  ',zmax/(sqrt(solsperyear)*delta),'times seasonal dry skin depth'
155 !    print *,'  ',zmax/(sqrt(solsperyear)*delta*stretch),'times seasonal filled skin depth'
156 !    print *,'  Initial ice depth=',zdepthT(k)
157 !    print *
158  enddo
159 ! call outputmoduleparameters
160 ! print *
161
162  ! open and name all output files
163!  open(unit=34,file='subout.'//ext,action='write',status='unknown')
164!  open(unit=36,file='depthF.'//ext,action='write',status='unknown')
165!  open(unit=37,file='depths.'//ext,action='write',status='unknown')
166
167 ! print *,'Equilibrating initial temperature'
168 ! do i=1,4
169 !    call icelayer_mars(0d0,nz,NP,thIn,rhoc,z,porosity,pfrost,Tb,zdepthF, &
170 !      &  zdepthE,porefill(1:nz,:),Tmean1,Tmean3,zdepthG, &
171 !      &  latitude,albedo,p0,ecc,omega,eps,icefrac,zdepthT,avrho1, &
172 !      &  avrho1prescribed)
173 ! enddo
174
175  !print *,'History begins here'
176 porefill(1:nz,1:NP) =  porefill_in(1:nz,1:NP)
177  zdepthT(1:NP) = ssi_depth_in
178  do
179  !print *,'Zt0=  ',ZdepthT
180     call icelayer_mars(timestep,nz,NP,thIn,rhoc,z,porosity,pfrost,Tb,zdepthF, &
181          & zdepthE,porefill,Tmean1,Tmean3,zdepthG, &
182          & albedo,p0,icefrac,zdepthT,avrho1, &
183          & avrho1prescribed)
184     icetime = icetime+timestep
185   !  print *,'T_after=  ',Tb(:)
186 !    print *,'z=  ',z(:)
187 !    print *,'Zt=  ',ZdepthT
188     ssi_depth=ZdepthT(1)
189    ! if (abs(mod(icetime/100.,1.d0))<1.e-3) then ! output every 1000 years
190    !    do k=1,NP
191         !write(36,*) icetime,latitude(k),zdepthF(k),porefill(1:nz,k)
192           ! compact output format
193 !          write(36,'(f10.0,2x,f7.3,1x,f11.5,1x)',advance='no') &
194 !               & icetime,latitude(k),zdepthF(k)
195 !          call compactoutput(36,porefill(:,k),nz)
196 !          write(37,501) icetime,latitude(k),zdepthT(k), &
197 !               & Tmean1(k),Tmean3(k),zdepthG(k),avrho1(k)
198  !      enddo
199  !   endif
200!     print *,icetime
201     if (icetime>=tlast) exit
202  enddo
203
204 ! close(34)
205 ! close(36); close(37)
206
207!501 format (f10.0,2x,f7.3,2x,f10.4,2(2x,f6.2),2x,f9.3,2x,g11.4)
208
209end subroutine dyn_ss_ice_m
210
211!=======================================================================
212
213subroutine soilthprop(porosity,fill,rhocobs,tiobs,layertype, &
214     &     newrhoc,newti,icefrac)
215!***********************************************************************
216! soilthprop: assign thermal properties of icy soil or dirty ice
217!
218!     porositiy = void space / total volume
219!     rhof = density of free ice in space not occupied by regolith [kg/m^3]
220!     fill = rhof/icedensity <=1 (only relevant for layertype 1)
221!     rhocobs = heat capacity per volume of dry regolith [J/m^3]
222!     tiobs = thermal inertia of dry regolith [SI-units]
223!     layertype: 1=interstitial ice, 2=pure ice or ice with dirt
224!                3=pure ice, 4=ice-cemented soil, 5=custom values
225!     icefrac: fraction of ice in icelayer (only relevant for layertype 2)
226!     output are newti and newrhoc
227!***********************************************************************
228
229! DEPENDENCIES
230! ------------
231use stoppage, only: stop_clean
232
233! DECLARATION
234! -----------
235implicit none
236  integer, intent(IN) :: layertype
237  real(8), intent(IN) :: porosity, fill, rhocobs, tiobs
238  real(8), intent(OUT) :: newti, newrhoc
239  real(8), intent(IN) :: icefrac
240  real(8) kobs, cice, icedensity, kice
241  !parameter (cice=2000.d0, icedensity=926.d0, kice=2.4d0) ! unaffected by scaling
242  parameter (cice=1540.d0, icedensity=927.d0, kice=3.2d0) ! at 198 Kelvin
243  real(8) fA, ki0, ki, k
244  real(8), parameter :: kw=3. ! Mellon et al., JGR 102, 19357 (1997)
245
246  kobs = tiobs**2/rhocobs
247  ! k, rhoc, and ti are defined in between grid points
248  ! rhof and T are defined on grid points
249
250  newrhoc = -9999.
251  newti  = -9999.
252
253  select case (layertype)
254  case (1) ! interstitial ice
255     newrhoc = rhocobs + porosity*fill*icedensity*cice
256     if (fill>0.) then
257        !--linear addition (option A)
258        k = porosity*fill*kice + kobs
259        !--Mellon et al. 1997 (option B)
260        ki0 = porosity/(1/kobs-(1-porosity)/kw)
261        fA = sqrt(fill)
262        ki = (1-fA)*ki0 + fA*kice
263        !k = kw*ki/((1-porosity)*ki+porosity*kw)
264     else
265        k = kobs
266     endif
267     newti = sqrt(newrhoc*k)
268
269  case (2)  ! massive ice (pure or dirty ice)
270     newrhoc = rhocobs*(1.-icefrac)/(1.-porosity) + icefrac*icedensity*cice
271     k = icefrac*kice + (1.-icefrac)*kw
272     newti = sqrt(newrhoc*k)
273
274  case (3)  ! all ice, special case of layertype 2, which doesn't use porosity
275     newrhoc = icedensity*cice
276     k = kice
277     newti = sqrt(newrhoc*k)
278
279  case (4)  ! pores completely filled with ice, special case of layertype 1
280     newrhoc = rhocobs + porosity*icedensity*cice
281     k = porosity*kice + kobs ! option A, end-member case of type 1, option A
282     !k = kw*kice/((1-porosity)*kice+porosity*kw) ! option B, harmonic average
283     newti = sqrt(newrhoc*k)
284
285  case (5)  ! custom values
286     ! values from Mellon et al. (2004) for ice-cemented soil
287     newrhoc = 2018.*1040.
288     k = 2.5
289     newti = sqrt(newrhoc*k)
290
291  case default
292     call stop_clean(__FILE__,__LINE__,'invalid layer type',1)
293
294  end select
295
296end subroutine soilthprop
297
298
299!=======================================================================
300
301      real*8 function frostpoint(p)
302!     inverse of psv
303!     input is partial pressure [Pascal]
304!     output is temperature [Kelvin]
305      ! DECLARATION
306! -----------
307implicit none
308      real*8 p
309
310!-----inverse of parametrization 1
311!      real*8 DHmelt,DHvap,DHsub,R,pt,Tt
312!      parameter (DHmelt=6008.,DHvap=45050.)
313!      parameter (DHsub=DHmelt+DHvap)
314!      parameter (R=8.314,pt=6.11e2,Tt=273.16)
315!      frostpoint = 1./(1./Tt-R/DHsub*log(p/pt))
316
317!-----inverse of parametrization 2
318!     inverse of eq. (2) in Murphy & Koop (2005)
319      real*8 A,B
320      parameter (A=-6143.7, B=28.9074)
321      frostpoint = A / (log(p) - B)
322
323!-----approximate inverse of parametrization 3
324!     eq. (8) in Murphy & Koop (2005)
325!      frostpoint = (1.814625*log(p) + 6190.134)/(29.120 - log(p))
326
327      end function frostpoint
328
329
330!=======================================================================
331!=========================== fast_subs_mars ============================
332!=======================================================================
333
334
335!*****************************************************
336! Subroutines for fast method
337! written by Norbert Schorghofer 2007-2011
338!*****************************************************
339
340
341subroutine icelayer_mars(bigstep,nz,NP,thIn,rhoc,z,porosity,pfrost, &
342     & Tb,zdepthF,zdepthE,porefill,Tmean1,Tmean3,zdepthG, &
343     & albedo,p0,icefrac,zdepthT,avrho1, &
344     & avrho1prescribed)
345!*************************************************************************
346! bigstep = time step [Earth years]
347! latitude  [degree]
348!*************************************************************************
349  !use ice_table,      only: rho_ice
350  !use omp_lib
351  ! DECLARATION
352! -----------
353implicit none
354  integer, intent(IN) :: nz, NP
355  real(8), intent(IN) :: bigstep
356  real(8), intent(IN) :: thIn(NP), rhoc(NP), z(NMAX), porosity, pfrost(NP)
357  real(8), intent(INOUT) :: porefill(nz,NP), zdepthF(NP), zdepthT(NP)
358  real(8), intent(INOUT) :: Tb(nz)
359  real(8), intent(OUT), dimension(NP) :: zdepthE, Tmean1, Tmean3, zdepthG
360  real(8), intent(IN), dimension(NP) ::  albedo, p0
361  real(8), intent(IN) ::  icefrac
362  real(8), intent(OUT) :: avrho1(NP)
363  real(8), intent(IN), optional :: avrho1prescribed(NP)  ! <0 means absent
364
365  integer k, typeF, typeG, typeT, j, jump, typeP
366  real(8) fracIR, fracDust, ti(NMAX), rhocv(NMAX)
367  real(8) Diff, ypp(nz), avdrho(NP), avdrhoP(NP), B, deltaz
368  real(8), SAVE :: avdrho_old(100), zdepth_old(100)  ! NP<=100
369  logical mode2
370
371  avdrho_old = 1. ! initialization
372
373  !$omp parallel &
374  !$omp    private(Diff,fracIR,fracDust,B,typeT,j,ti,rhocv,typeF,jump,typeG)
375  !$omp do
376  do k=1,NP   ! big loop
377
378     Diff = 4e-4*600./p0(k)
379     fracIR = 0.04*p0(k)/600.; fracDust = 0.02*p0(k)/600.
380     B = Diff*bigstep*86400.*365.24/(porosity*927.)
381     !B = Diff*bigstep*86400.*365.24/(porosity*rho_ice(Tb(),'H2O'))
382
383     typeT = -9
384     if (zdepthT(k)>=0. .and. zdepthT(k)<z(nz)) then
385        do j=1,nz
386           if (z(j)>zdepthT(k)) then ! ice
387              typeT = j  ! first point in ice
388              exit
389           endif
390        enddo
391     endif
392
393   !  call assignthermalproperties(nz,thIn(k),rhoc(k),ti,rhocv, &
394   !       & typeT,icefrac,porosity,porefill(:,k))
395
396     !----run thermal model
397     call ajsub_mars(typeT, albedo(k), pfrost(k), nz, z, &
398          &     ti, rhocv, fracIR, fracDust, p0(k), &
399          &     avdrho(k), avdrhoP(k), avrho1(k), Tb(k), zdepthE(k), typeF, &
400          &     zdepthF(k), ypp, porefill(:,k), Tmean1(k), Tmean3(k), B, &
401          &     typeG, zdepthG(k), avrho1prescribed(k))
402
403     if (typeF*zdepthF(k)<0.) stop 'error in icelayer_mars'
404     ! diagnose
405     if (zdepthT(k)>=0.) then
406        jump = 0
407        do j=1,nz
408           if (zdepth_old(k)<z(j) .and. zdepthT(k)>z(j)) jump = jump+1
409        enddo
410     else
411        jump = -9
412     endif
413     if (zdepthT(k)>=0. .and. avdrho(k)*avdrho_old(k)<0.) then
414        write(34,*) '# zdepth arrested'
415        if (jump>1) write(34,*) '# previous step was too large',jump
416     endif
417!     write(34,'(f8.3,1x,f6.2,1x,f11.5,2(1x,g11.4),1x,i3)') &
418!          &        bigstep,latitude(k),zdepthT(k),avdrho(k),avdrhoP(k),jump
419     zdepth_old(k) = zdepthT(k)
420     avdrho_old(k) = avdrho(k)
421
422!----mode 2 growth
423     typeP = -9;  mode2 = .FALSE.
424     do j=1,nz
425        if (porefill(j,k)>0.) then
426           typeP = j  ! first point with ice
427           exit
428        endif
429     enddo
430     if (typeT>0 .and. typeP>2 .and. zdepthE(k)>0.) then
431        if (porefill(typeP,k)>=1. .and. porefill(typeP-1,k)==0. .and. &
432             & zdepthE(k)<z(typeP) .and. &
433             & z(typeP)-zdepthE(k)>2*(z(typeP)-z(typeP-1))) then  ! trick that avoids oscillations
434           deltaz = -avdrhoP(k)/z(typeP)*18./8314.*B  ! conservation of mass
435           if (deltaz>z(typeP)-z(typeP-1)) then  ! also implies avdrhoP<0.
436              mode2 = .TRUE.
437           endif
438        endif
439     endif
440
441     !call icechanges_poreonly(nz,z,typeF,typeG,avdrhoP(k),ypp,B,porefill(:,k))
442     call icechanges(nz,z(:),typeF,avdrho(k),avdrhoP(k),ypp(:), &
443          & Diff,porosity,icefrac,bigstep,zdepthT(k),porefill(:,k),typeG)
444     if (typeP>2) then
445     if (mode2 .and. porefill(typeP,k)>=1. .and. porefill(typeP-1,k)==0.) then  ! nothing changed
446        porefill(typeP-1,k)=1.  ! paste a layer
447   !     write(34,*) '# mode 2 growth occurred',typeP,typeF,typeT
448     endif
449     endif
450
451  enddo  ! end of big loop
452  !$omp end do
453  !$omp end parallel
454end subroutine icelayer_mars
455
456
457
458subroutine ajsub_mars(typeT, albedo0, pfrost, nz, z, ti, rhocv, &
459     &     fracIR, fracDust, patm, avdrho, avdrhoP, avrho1, &
460     &     Tb, zdepthE, typeF, zdepthF, ypp, porefill, Tmean1, Tmean3, &
461     &     B, typeG, zdepthG, avrho1prescribed)
462!***********************************************************************
463!  A 1D thermal model that returns various averaged quantities
464!
465!  Tb<0 initializes temperatures
466!  Tb>0 initializes temperatures with Tb
467!***********************************************************************
468  use orbit, only: yr_len_sol, sol_len_s
469  ! DECLARATION
470! -----------
471implicit none
472  integer, intent(IN) :: nz, typeT
473!  real(8), intent(IN) :: latitude  ! in radians
474  real(8), intent(IN) :: albedo0, pfrost, z(NMAX)
475  real(8), intent(IN) :: ti(NMAX), rhocv(NMAX), fracIR, fracDust, patm
476  real(8), intent(IN) ::  porefill(nz)
477  real(8), intent(OUT) :: avdrho, avdrhoP  ! difference in annual mean vapor density
478  real(8), intent(OUT) :: avrho1  ! mean vapor density on surface
479  real(8), intent(INOUT) :: Tmean1
480  real(8), intent(INOUT) :: Tb(nz)
481  integer, intent(OUT) :: typeF  ! index of depth below which filling occurs
482  real(8), intent(OUT) :: zdepthE, zdepthF
483  real(8), intent(OUT) :: ypp(nz) ! (d rho/dt)/Diff
484  real(8), intent(OUT) :: Tmean3, zdepthG
485  real(8), intent(IN) :: B  ! just passing through
486  integer, intent(OUT) :: typeG
487  real(8), intent(IN), optional :: avrho1prescribed
488  real(8), parameter :: a = 1.52366 ! Mars semimajor axis in A.U.
489  integer nsteps, n, i, nm, typeP
490  real(8) tmax, time, Qn, Qnp1, tdays
491  real(8) marsR, marsLs, marsDec, HA
492  real(8) Tsurf, Tco2frost, albedo, Fsurf, m, dE, emiss, T(NMAX)
493  real(8) Told(nz), Fsurfold, Tsurfold, Tmean0, avrho2
494  real(8) rhosatav0, rhosatav(nz), rlow
495
496  tmax = EQUILTIME*yr_len_sol
497  nsteps = int(tmax/dt)     ! calculate total number of timesteps
498
499  Tco2frost = tfrostco2(patm)
500
501  !if (Tb<=0.) then  ! initialize
502     !Tmean0 = 210.15       ! black-body temperature of planet
503     !Tmean0 = (589.*(1.-albedo0)*cos(latitude)/(pi*emiss0*sigSB))**0.25 ! better estimate
504     !Tmean0 = Tmean0-5.
505     !write(34,*) '# initialized with temperature estimate at',latitude/d2r,'of',Tmean0,'K'
506     !T(1:nz) = Tmean0
507  !else
508     T(1:nz) = Tb
509     ! not so good when Fgeotherm is on
510  !endif
511
512  albedo = albedo0
513  emiss = emiss0
514  do i=1,nz
515     if (T(i)<Tco2frost) T(i)=Tco2frost
516  enddo
517  Tsurf = T(1)
518  m=0.; Fsurf=0.
519
520  nm=0
521  avrho1=0.; avrho2=0.
522  Tmean1=0.; Tmean3=0.
523  rhosatav0 = 0.
524  rhosatav(:) = 0.
525
526  time=0.
527  !call generalorbit(0.d0,a,ecc,omega,eps,marsLs,marsDec,marsR)
528  !HA = 2.*pi*time            ! hour angle
529!  Qn=flux(marsR,marsDec,latitude,HA,albedo,fracir,fracdust,0.d0,0.d0)
530  !Qn = flux_mars77(marsR,marsDec,latitude,HA,albedo,fracir,fracdust)
531  !----loop over time steps
532  do n=0,nsteps-1
533     time = (n+1)*dt         !   time at n+1
534   !  tdays = time*(sol_len_s/earthDay) ! parenthesis may improve roundoff
535   !  call generalorbit(tdays,a,ecc,omega,eps,marsLs,marsDec,marsR)
536   !  HA = 2.*pi*mod(time,1.d0)  ! hour angle
537!     Qnp1=flux(marsR,marsDec,latitude,HA,albedo,fracir,fracdust,0.d0,0.d0)
538     !Qnp1 = flux_mars77(marsR,marsDec,latitude,HA,albedo,fracir,fracdust)
539
540     Tsurfold = Tsurf
541     Fsurfold = Fsurf
542     Told(1:nz) = T(1:nz)
543     if (m<=0. .or. Tsurf>Tco2frost) then
544       ! call conductionQ(nz,z,dt*sol_len_s,Qn,Qnp1,T,ti,rhocv,emiss, &
545       !      &           Tsurf,Fgeotherm,Fsurf)
546     endif
547     if (Tsurf<Tco2frost .or. m>0.) then ! CO2 condensation
548        T(1:nz) = Told(1:nz)
549        !call conductionT(nz,z,dt*sol_len_s,T,Tsurfold,Tco2frost,ti, &
550             !&              rhocv,Fgeotherm,Fsurf)
551        Tsurf = Tco2frost
552    !    dE = (- Qn - Qnp1 + Fsurfold + Fsurf + &
553    !         &      emiss*sigSB*(Tsurfold**4+Tsurf**4)/2.
554        m = m + dt*sol_len_s*dE/Lco2frost
555     endif
556     if (Tsurf>Tco2frost .or. m<=0.) then
557        albedo = albedo0
558        emiss = emiss0
559     else
560        albedo = co2albedo
561        emiss = co2emiss
562     endif
563     !Qn=Qnp1
564
565     if (time>=tmax-yr_len_sol) then
566        Tmean1 = Tmean1 + Tsurf
567        Tmean3 = Tmean3 + T(nz)
568        avrho1 = avrho1 + min(psv(Tsurf),pfrost)/Tsurf
569        rhosatav0 = rhosatav0 + psv(Tsurf)/Tsurf
570        do i=1,nz
571           rhosatav(i) = rhosatav(i) + psv(T(i))/T(i)
572        enddo
573        nm = nm+1
574     endif
575
576  enddo  ! end of time loop
577
578  avrho1 = avrho1/nm
579  if (present(avrho1prescribed)) then
580     if (avrho1prescribed>=0.) avrho1=avrho1prescribed
581  endif
582  Tmean1 = Tmean1/nm; Tmean3 = Tmean3/nm
583  rhosatav0 = rhosatav0/nm; rhosatav(:)=rhosatav(:)/nm
584  if (typeT>0) then
585     avrho2 = rhosatav(typeT)
586  else
587     avrho2 = rhosatav(nz) ! no ice
588  endif
589  avdrho = avrho2-avrho1
590  typeP = -9
591  do i=1,nz
592     if (porefill(i)>0.) then
593        typeP = i  ! first point with ice
594        exit
595     endif
596  enddo
597  if (typeP>0) then
598     avdrhoP = rhosatav(typeP) - avrho1
599  else
600     avdrhoP = -9999.
601  end if
602
603  zdepthE = equildepth(nz, z, rhosatav, rhosatav0, avrho1)
604
605  if (Fgeotherm>0.) then
606     !Tb = Tmean1
607     typeG = 1   ! will be overwritten by depths_avmeth
608     rlow = 2*rhosatav(nz)-rhosatav(nz-1)
609  else
610     !Tb = T(nz)
611     typeG = -9
612     rlow = rhosatav(nz-1)
613  endif
614  call depths_avmeth(typeT, nz, z, rhosatav(:), rhosatav0, rlow, avrho1,  &
615       & porefill(:), typeF, zdepthF, B, ypp(:), typeG, zdepthG)
616
617end subroutine ajsub_mars
618
619
620
621real*8 function tfrostco2(p)
622!     the inverse of function psvco2
623!     input is pressure in Pascal, output is temperature in Kelvin
624! DECLARATION
625! -----------
626implicit none
627real*8 p
628tfrostco2 = 3182.48/(23.3494+log(100./p))
629end function
630
631!=======================================================================
632
633      real*8 function psv(T)
634!     saturation vapor pressure of H2O ice [Pascal]
635!     input is temperature [Kelvin]
636      ! DECLARATION
637! -----------
638implicit none
639      real*8 T
640
641!-----parametrization 1
642!      real*8 DHmelt,DHvap,DHsub,R,pt,Tt,C
643!      parameter (DHmelt=6008.,DHvap=45050.)
644!      parameter (DHsub=DHmelt+DHvap) ! sublimation enthalpy [J/mol]
645!      parameter (R=8.314,pt=6.11e2,Tt=273.16)
646!      C = (DHsub/R)*(1./T - 1./Tt)
647!      psv = pt*exp(-C)
648
649!-----parametrization 2
650!     eq. (2) in Murphy & Koop, Q. J. R. Meteor. Soc. 131, 1539 (2005)
651!     differs from parametrization 1 by only 0.1%
652      real*8 A,B
653      parameter (A=-6143.7, B=28.9074)
654      psv = exp(A/T+B)  ! Clapeyron
655
656!-----parametrization 3
657!     eq. (7) in Murphy & Koop, Q. J. R. Meteor. Soc. 131, 1539 (2005)
658!     psv = exp(9.550426 - 5723.265/T + 3.53068*log(T) - 0.00728332*T)
659
660      end function psv
661
662
663!=======================================================================
664!=========================== fast_subs_univ ============================
665!=======================================================================
666
667
668!*****************************************************
669! Commonly used subroutines for fast method
670! written by Norbert Schorghofer 2007-2011
671!*****************************************************
672
673pure function zint(y1,y2,z1,z2)
674  ! interpolate between two values, y1*y2<0
675  ! DECLARATION
676! -----------
677implicit none
678  real(8), intent(IN) :: y1,y2,z1,z2
679  real(8) zint
680  zint = (y1*z2 - y2*z1)/(y1-y2)
681end function zint
682
683
684
685pure function equildepth(nz, z, rhosatav, rhosatav0, avrho1)
686!***********************************************************************
687!  returns equilibrium depth for given ice content
688!  this is not the true (final) equilibrium depth
689!***********************************************************************
690  ! DECLARATION
691! -----------
692implicit none
693  integer, intent(IN) :: nz
694  real(8), intent(IN) :: z(nz), rhosatav(nz)
695  real(8), intent(IN) :: rhosatav0, avrho1
696  integer i, typeE
697  real(8) equildepth  ! =zdepthE
698  !real(8), external :: zint  ! defined in allinterfaces.mod
699
700  typeE = -9; equildepth = -9999.
701  do i=1,nz
702     if (rhosatav(i) <= avrho1) then
703        typeE=i
704        exit
705     endif
706  enddo
707  if (typeE>1) then ! interpolate
708     equildepth=zint(avrho1-rhosatav(typeE-1),avrho1-rhosatav(typeE),z(typeE-1),z(typeE))
709  end if
710  if (typeE==1) equildepth=zint(avrho1-rhosatav0,avrho1-rhosatav(1),0.d0,z(1))
711  if (equildepth>z(nz)) equildepth=-9999.   ! happens very rarely
712end function equildepth
713
714
715
716subroutine depths_avmeth(typeT, nz, z, rhosatav, rhosatav0, rlow, avrho1,  &
717     & porefill, typeF, zdepthF, B, ypp, typeG, zdepthG)
718!***********************************************************************
719!  returns interface depth and ypp
720!  also returns lower geothermally limited boundary, if applicable
721!
722!  this is a crucial subroutine
723!
724!  B = Diff*bigstep/(porosity*icedensity)  [SI units]
725!***********************************************************************
726  use maths,     only: deriv2_simple, deriv1_onesided, deriv1, colint
727  ! DECLARATION
728! -----------
729implicit none
730  integer, intent(IN) :: nz, typeT
731  real(8), intent(IN), dimension(nz) :: z, rhosatav, porefill
732  real(8), intent(IN) :: rhosatav0, rlow, avrho1
733  integer, intent(OUT) :: typeF  ! index of depth below which filling occurs
734  real(8), intent(INOUT) :: zdepthF
735  real(8), intent(IN) :: B
736  real(8), intent(OUT) :: ypp(nz), zdepthG
737  integer, intent(INOUT) :: typeG  ! positive on input when Fgeotherm>0
738
739  integer i, typeP, nlb, newtypeG
740  real(8) eta(nz), Jpump1, help(nz), yp(nz), zdepthFold, ap_one, ap(nz)
741  real(8) leak, cumfill, cumfillabove
742
743  if (typeT<0) then
744     nlb = nz
745     do i=1,nz
746        eta(i) = constriction(porefill(i))
747     enddo
748  else
749     !nlb = typeT-1
750     nlb = typeT ! changed 2010-09-29
751     do i=1,typeT-1
752        eta(i) = constriction(porefill(i))
753     enddo
754     do i=typeT,nz
755        eta(i)=0.
756     enddo
757  end if
758
759!-fill depth
760  zdepthFold = zdepthF
761  typeF = -9;  zdepthF = -9999.
762  call deriv1(z,nz,rhosatav,rhosatav0,rlow,yp)  ! yp also used below
763  do i=1,nlb
764     Jpump1 = (rhosatav(i)-avrho1)/z(i)  ! <0 when stable
765     ! yp is always <0
766     help(i) = Jpump1 - eta(i)*yp(i)
767     leak = porefill(i)/B*(z(i)-zdepthFold)/(18./8314.)
768     !help(i) = help(i)-leak   ! optional
769     if (help(i) <= 0.) then
770        typeF=i
771        !print *,'#',typeF,Jpump1,eta(typeF)*yp(typeF),leak
772        exit
773     endif
774  enddo
775  if (typeF>1) zdepthF = zint(help(typeF-1),help(typeF),z(typeF-1),z(typeF))
776  if (typeF==1) zdepthF=z(1)
777
778
779!-depth to shallowest perennial ice
780  typeP = -9
781  do i=1,nz
782     if (porefill(i)>0.) then
783        typeP = i  ! first point with ice
784        exit
785     endif
786  enddo
787
788!-calculate ypp
789  !call deriv1(z,nz,rhosatav,rhosatav0,rlow,yp)
790  call deriv1(z,nz,eta(:),1.d0,eta(nz-1),ap)
791  if (typeP>0 .and. typeP<nz-2) then
792     call deriv1_onesided(typeP,z,nz,eta(:),ap_one)
793     ! print *,typeP,ap(typeP),ap_one
794     ap(typeP)=ap_one
795  endif
796  call deriv2_simple(z,nz,rhosatav(1:nz),rhosatav0,rlow,ypp(:))
797  !call deriv2_full(z,nz,eta(:),rhosatav(:),1.d0,rhosatav0,rhosatav(nz-1),ypp(:))
798  !write(40,*) rhosatav
799  !write(41,*) yp
800  !write(42,*) ypp
801  ypp(:) = ap(:)*yp(1:)+eta(:)*ypp(:)
802  !write(43,*) ypp
803  !write(44,*) eta(1:nz)
804  !write(45,*) (rhosatav(:)-avrho1)/z(:)
805  ypp(:) = ypp(:)*18./8314.
806  ! ypp values beyond nlb should not be used
807
808!-geothermal stuff
809  if (typeT>0) typeG=-9
810  if (typeG<0) zdepthG=-9999.
811  if (typeG>0 .and. typeT<0) then
812     typeG=-9
813     do i=2,nz
814        if (yp(i)>0.) then  ! first point with reversed flux
815           typeG=i
816           zdepthG=zint(yp(i-1),yp(i),z(i-1),z(i))
817           !zdepthG=z(typeG)
818           exit
819        endif
820     enddo
821  else
822     typeG = -9
823  endif
824  if (typeG>0 .and. typeT<0) then
825     call colint(porefill(:)/eta(:),z,nz,typeG-1,nz,cumfillabove)
826     newtypeG = -9
827     do i=typeG,nz
828        if (minval(eta(i:nz))<=0.) cycle  ! eta=0 means completely full
829        call colint(porefill(:)/eta(:),z,nz,i,nz,cumfill)
830        if (cumfill<yp(i)*18./8314.*B) then  ! usually executes on i=typeG
831           if (i>typeG) then
832              write(34,*) '# adjustment to geotherm depth by',i-typeG
833              zdepthG = zint(yp(i-1)*18./8314.*B-cumfillabove, &  ! no good
834                   &        yp(i)*18./8314.*B-cumfill,z(i-1),z(i))
835              if (zdepthG>z(i) .or. zdepthG<z(i-1)) write(34,*) &
836                   & '# WARNING: zdepthG interpolation failed',i,z(i-1),zdepthG,z(i)
837              newtypeG=i
838           endif
839           ! otherwise leave zdepthG untouched
840           exit
841        endif
842        cumfillabove = cumfill
843     enddo
844     if (newtypeG>0) typeG=newtypeG
845  end if
846  ! if typeG>0, then all ice at and below typeG should be erased
847end subroutine depths_avmeth
848
849
850
851pure subroutine icechanges(nz,z,typeF,avdrho,avdrhoP,ypp, &
852     & Diff,porosity,icefrac,bigstep,zdepthT,porefill,typeG)
853!***********************************************************
854! advances ice table, advances interface, grows pore ice
855!
856! a crucial subroutine
857!***********************************************************
858  use maths,     only: colint
859  !use ice_table, only: rho_ice
860  ! DECLARATION
861! -----------
862implicit none
863  integer, intent(IN) :: nz, typeF, typeG
864  real(8), intent(IN) :: z(nz), ypp(nz), avdrho, avdrhoP
865  real(8), intent(IN) :: Diff, porosity, icefrac, bigstep
866  real(8), intent(INOUT) :: zdepthT, porefill(nz)
867  integer j, erase, newtypeP, ub, typeP, typeT
868  real(8) B, beta, integ
869
870  B = Diff*bigstep*86400.*365.24/(porosity*927.)
871  !B = Diff*bigstep*86400.*365.24/(porosity*rho_ice(T,'H2O'))
872
873  ! advance ice table, avdrho>0 is retreat
874  if (zdepthT>=0. .and. avdrho>0.) then
875     typeP=-9999; typeT=-9999
876     do j=1,nz
877        if (z(j)>zdepthT) then
878           typeT=j
879           exit
880        endif
881     enddo
882     do j=1,nz
883        if (porefill(j)>0.) then
884           typeP=j
885           exit
886        endif
887     enddo
888     if (typeP==typeT) then   ! new 2011-09-01
889        beta = (1-icefrac)/(1-porosity)/icefrac
890        beta = Diff*bigstep*beta*86400*365.24/927.
891        !beta = Diff*bigstep*beta*86400*365.24/rho_ice(T,'H2O')
892        zdepthT = sqrt(2*beta*avdrho*18./8314. + zdepthT**2)
893     endif
894  endif
895  if (zdepthT>z(nz)) zdepthT=-9999.
896
897  ! advance interface, avdrhoP>0 is loss from zdepthP
898  if (avdrhoP>0.) then
899     erase=0
900     do j=1,nz
901        if (typeF>0 .and. j>=typeF) exit ! don't retreat beyond typeF
902        if (zdepthT>=0. .and. z(j)>zdepthT) exit
903        call colint(porefill(1:nz)*z(1:nz),z(1:nz),nz,1,j,integ)
904        erase = j
905        if (integ>B*avdrhoP*18./8314.) exit
906     end do
907     if (erase>0) porefill(1:erase)=0.
908  endif
909
910  ! new depth
911  newtypeP = -9
912  do j=1,nz
913     if (zdepthT>=0. .and. z(j)>zdepthT) exit
914     if (porefill(j)>0.) then
915        newtypeP = j  ! first point with pore ice
916        exit
917     endif
918  enddo
919
920  ! diffusive filling
921  ub = typeF
922  if (newtypeP>0 .and. typeF>0 .and. newtypeP<ub) ub=newtypeP
923  if (ub>0) then
924     do j=ub,nz
925        porefill(j) = porefill(j) + B*ypp(j)
926        if (porefill(j)<0.) porefill(j)=0.
927        if (porefill(j)>1.) porefill(j)=1.
928        if (zdepthT>=0. .and. z(j)>zdepthT) exit
929     enddo
930  end if
931
932  ! below icetable
933  if (zdepthT>=0.) then
934     do j=1,nz
935        if (z(j)>zdepthT) porefill(j) = icefrac/porosity
936     enddo
937  else
938     ! geothermal lower boundary
939     if (typeG>0) porefill(typeG:nz)=0.
940  end if
941end subroutine icechanges
942
943!=======================================================================
944FUNCTION constriction(porefill) RESULT(eta)
945!-----------------------------------------------------------------------
946! NAME
947!     constriction
948!
949! DESCRIPTION
950!     Compute the constriction of vapor flux by pore ice.
951!
952! AUTHORS & DATE
953!     L. Lange
954!     JB Clement, 2023-2025
955!
956! NOTES
957!
958!-----------------------------------------------------------------------
959
960! DECLARATION
961! -----------
962implicit none
963
964! ARGUMENTS
965! ---------
966real, intent(in) :: porefill ! pore filling fraction
967
968! LOCAL VARIABLES
969! ---------------
970real :: eta ! constriction
971
972! CODE
973! ----
974if (porefill <= 0.) then
975    eta = 1.
976else if (0 < porefill .and. porefill < 1.) then
977    eta = (1-porefill)**2 ! Hudson et al., JGR, 2009
978else
979    eta = 0.
980endif
981
982END FUNCTION constriction
983!=======================================================================
984
985END MODULE subsurface_ice
Note: See TracBrowser for help on using the repository browser.