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