source: trunk/LMDZ.COMMON/libf/evolution/NS_dyn_ss_ice_m.F90 @ 3527

Last change on this file since 3527 was 3527, checked in by jbclement, 26 hours ago

PEM:

  • Removing redundant Norbert Schorghofer's subroutines/parameters (follow-up of r3526);
  • Making all Norbert Schorghofer's subroutines with modern explicit interface via modules;
  • Cleaning of "glaciers_mod.F90";
  • Optimization for the computation of ice density according to temperature by using a function.

JBC

File size: 9.8 KB
Line 
1MODULE dyn_ss_ice_m_mod
2
3implicit none
4
5CONTAINS
6
7!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
8!!!
9!!! Purpose:  Retreat and growth of subsurface ice on Mars
10!!!           orbital elements remain constant
11!!!
12!!!
13!!! Author: EV, updated NS MSIM dynamical program for the PEM
14!!!
15!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
16
17SUBROUTINE dyn_ss_ice_m(ssi_depth_in,T1,Tb,nz,thIn,p0,pfrost,porefill_in,porefill,ssi_depth)
18
19!***********************************************************************
20! Retreat and growth of subsurface ice on Mars
21! orbital elements remain constant
22!***********************************************************************
23  use constants_marspem_mod, only: sec_per_sol
24  use fast_subs_mars, only: psv, icelayer_mars, NMAX
25#ifndef CPP_STD
26    use comcstfi_h,   only: pi
27#else
28    use comcstfi_mod, only: pi
29#endif
30  implicit none
31  integer, parameter :: NP=1   ! # of sites
32  integer nz, i, k, iloop
33  real(8)  zmax, delta, z(NMAX), icetime, porosity, icefrac
34  real(8), dimension(NP) ::  albedo, thIn, rhoc
35  real(8), dimension(NP) :: pfrost, p0
36  real(8) newti, stretch, newrhoc, ecc, omega, eps, timestep
37  real(8)  ssi_depth_in, ssi_depth, T1
38  real(8), dimension(NP) :: zdepthF, zdepthE, zdepthT, zdepthG
39  real(8), dimension(NMAX,NP) :: porefill, porefill_in
40  real(8), dimension(nz) :: Tb
41  real(8), dimension(NP) :: Tmean1, Tmean3, avrho1
42  real(8) tmax, tlast, avrho1prescribed(NP), l1
43  real(8), external :: smartzfac
44
45  !if (iargc() /= 1) then
46  !   stop 'USAGE: icages ext'
47  !endif
48  !call getarg( 1, ext )
49
50  if (NP>100) stop 'subroutine icelayer_mars cannot handle this many sites'
51
52  ! parameters that never ever change
53  porosity = 0.4d0  ! porosity of till
54  !rhoc(:) = 1500.*800.  ! will be overwritten
55  icefrac = 0.98
56  tmax = 1
57  tlast = 0.
58  avrho1prescribed(:) = pfrost/T1  ! <0 means absent
59  albedo=0.23
60  !avrho1prescribed(:) = 0.16/200.  ! units are Pa/K
61
62  !open(unit=21,file='lats.'//ext,action='read',status='old',iostat=ierr)
63  !if (ierr /= 0) then
64  !   print *,'File lats.'//ext,'not found'
65  !   stop
66  !endif
67  do k=1,NP
68     !read(21,*) latitude(k),albedo(k),thIn(k),htopo(k)
69     ! empirical relation from Mellon & Jakosky
70     rhoc(k) = 800.*(150.+100.*sqrt(34.2+0.714*thIn(k)))
71  enddo
72  !close(21)
73
74  ! set eternal grid
75  zmax = 25.
76  !zfac = smartzfac(nz,zmax,6,0.032d0)
77  !call setgrid(nz,z,zmax,zfac)
78  l1=2.e-4
79  do iloop=0,nz - 1
80    z(iloop + 1) = l1*(1+iloop**2.9*(1-exp(-real(iloop)/20.)))
81  enddo
82 
83
84  !open(unit=30,file='z.'//ext,action='write',status='unknown')
85  !write(30,'(999(f8.5,1x))') z(1:nz)
86  !close(30)
87
88  !ecc = ecc_in;  eps = obl_in*d2r;  omega = Lp_in*d2r   ! today
89  ! total atmospheric pressure
90  !p0(:) = 600.
91  ! presently 520 Pa at zero elevation (Smith & Zuber, 1998)
92 ! do k=1,NP
93 !    p0(k)=520*exp(-htopo(k)/10800.)
94 ! enddo
95  timestep = 1 ! must be integer fraction of 1 ka
96  icetime = -tmax-timestep  ! earth years
97 
98  ! initializations
99  !Tb = -9999.
100  zdepthF(:) = -9999.
101
102  !zdepthT(1:NP) = -9999.   ! reset again below
103 ! zdepthT(1:NP) = 0.
104
105 ! print *,'RUNNING MARS_FAST'
106 ! print *,'Global model parameters:'
107 ! print *,'nz=',nz,' zfac=',zfac,'zmax=',zmax
108 ! print *,'porosity=',porosity
109 ! print *,'starting at time',icetime,'years'
110 ! print *,'time step=',timestep,'years'
111 ! print *,'eps=',eps/d2r,'ecc=',ecc,'omega=',omega/d2r
112 ! print *,'number of sites=',NP
113 ! print *,'Site specific parameters:'
114  do k=1,NP
115     if (NP>1) print *,'  Site ',k
116 !    print *,'  latitude (deg)',latitude(k),' rho*c (J/m^3/K)',rhoc(k),' thIn=',thIn(k)
117 !    print *,'  total pressure=',p0(k),'partial pressure=',pfrost(k)
118     delta = thIn(k)/rhoc(k)*sqrt(sec_per_sol/pi)
119 !    print *,'  skin depths (m)',delta,delta*sqrt(solsperyear)
120     call soilthprop(porosity,1.d0,rhoc(k),thIn(k),1,newrhoc,newti,icefrac)
121     stretch = (newti/thIn(k))*(rhoc(k)/newrhoc)
122     do i=1,nz
123        if (z(i)<delta) cycle
124 !       print *,'  ',i-1,' grid points within diurnal skin depth'
125        exit
126     enddo
127 !    print *,'  ',zmax/(sqrt(solsperyear)*delta),'times seasonal dry skin depth'
128 !    print *,'  ',zmax/(sqrt(solsperyear)*delta*stretch),'times seasonal filled skin depth'
129 !    print *,'  Initial ice depth=',zdepthT(k)
130 !    print *
131  enddo
132 ! call outputmoduleparameters
133 ! print *
134
135  ! open and name all output files
136!  open(unit=34,file='subout.'//ext,action='write',status='unknown')
137!  open(unit=36,file='depthF.'//ext,action='write',status='unknown')
138!  open(unit=37,file='depths.'//ext,action='write',status='unknown')
139
140 ! print *,'Equilibrating initial temperature'
141 ! do i=1,4
142 !    call icelayer_mars(0d0,nz,NP,thIn,rhoc,z,porosity,pfrost,Tb,zdepthF, &
143 !      &  zdepthE,porefill(1:nz,:),Tmean1,Tmean3,zdepthG, &
144 !      &  latitude,albedo,p0,ecc,omega,eps,icefrac,zdepthT,avrho1, &
145 !      &  avrho1prescribed)
146 ! enddo
147
148  !print *,'History begins here'
149 porefill(1:nz,1:NP) =  porefill_in(1:nz,1:NP)
150  zdepthT(1:NP) = ssi_depth_in
151  do
152  !print *,'Zt0=  ',ZdepthT
153     call icelayer_mars(timestep,nz,NP,thIn,rhoc,z,porosity,pfrost,Tb,zdepthF, &
154          & zdepthE,porefill(1:nz,:),Tmean1,Tmean3,zdepthG, &
155          & albedo,p0,icefrac,zdepthT,avrho1, &
156          & avrho1prescribed)
157     icetime = icetime+timestep
158   !  print *,'T_after=  ',Tb(:)
159 !    print *,'z=  ',z(:)
160 !    print *,'Zt=  ',ZdepthT
161     ssi_depth=ZdepthT(1)
162    ! if (abs(mod(icetime/100.,1.d0))<1.e-3) then ! output every 1000 years
163    !    do k=1,NP
164         !write(36,*) icetime,latitude(k),zdepthF(k),porefill(1:nz,k)
165           ! compact output format
166 !          write(36,'(f10.0,2x,f7.3,1x,f11.5,1x)',advance='no') &
167 !               & icetime,latitude(k),zdepthF(k)
168 !          call compactoutput(36,porefill(:,k),nz)
169 !          write(37,501) icetime,latitude(k),zdepthT(k), &
170 !               & Tmean1(k),Tmean3(k),zdepthG(k),avrho1(k)
171  !      enddo
172  !   endif
173!     print *,icetime
174     if (icetime>=tlast) exit
175  enddo
176
177 ! close(34)
178 ! close(36); close(37)
179
180!501 format (f10.0,2x,f7.3,2x,f10.4,2(2x,f6.2),2x,f9.3,2x,g11.4)
181 
182end subroutine dyn_ss_ice_m
183
184!=======================================================================
185
186subroutine soilthprop(porosity,fill,rhocobs,tiobs,layertype, &
187     &     newrhoc,newti,icefrac)
188!***********************************************************************
189! soilthprop: assign thermal properties of icy soil or dirty ice
190!
191!     porositiy = void space / total volume
192!     rhof = density of free ice in space not occupied by regolith [kg/m^3]
193!     fill = rhof/icedensity <=1 (only relevant for layertype 1)
194!     rhocobs = heat capacity per volume of dry regolith [J/m^3]
195!     tiobs = thermal inertia of dry regolith [SI-units]
196!     layertype: 1=interstitial ice, 2=pure ice or ice with dirt
197!                3=pure ice, 4=ice-cemented soil, 5=custom values
198!     icefrac: fraction of ice in icelayer (only relevant for layertype 2)
199!     output are newti and newrhoc
200!***********************************************************************
201  implicit none
202  integer, intent(IN) :: layertype
203  real(8), intent(IN) :: porosity, fill, rhocobs, tiobs
204  real(8), intent(OUT) :: newti, newrhoc
205  real(8), intent(IN) :: icefrac
206  real(8) kobs, cice, icedensity, kice
207  !parameter (cice=2000.d0, icedensity=926.d0, kice=2.4d0) ! unaffected by scaling
208  parameter (cice=1540.d0, icedensity=927.d0, kice=3.2d0) ! at 198 Kelvin
209  real(8) fA, ki0, ki, k
210  real(8), parameter :: kw=3. ! Mellon et al., JGR 102, 19357 (1997)
211
212  kobs = tiobs**2/rhocobs
213  ! k, rhoc, and ti are defined in between grid points
214  ! rhof and T are defined on grid points
215
216  newrhoc = -9999.
217  newti  = -9999.
218
219  select case (layertype)
220  case (1) ! interstitial ice
221     newrhoc = rhocobs + porosity*fill*icedensity*cice
222     if (fill>0.) then
223        !--linear addition (option A)
224        k = porosity*fill*kice + kobs
225        !--Mellon et al. 1997 (option B)
226        ki0 = porosity/(1/kobs-(1-porosity)/kw)
227        fA = sqrt(fill)
228        ki = (1-fA)*ki0 + fA*kice
229        !k = kw*ki/((1-porosity)*ki+porosity*kw)
230     else
231        k = kobs
232     endif
233     newti = sqrt(newrhoc*k)
234     
235  case (2)  ! massive ice (pure or dirty ice)
236     newrhoc = rhocobs*(1.-icefrac)/(1.-porosity) + icefrac*icedensity*cice
237     k = icefrac*kice + (1.-icefrac)*kw
238     newti = sqrt(newrhoc*k)
239 
240  case (3)  ! all ice, special case of layertype 2, which doesn't use porosity
241     newrhoc = icedensity*cice
242     k = kice
243     newti = sqrt(newrhoc*k)
244 
245  case (4)  ! pores completely filled with ice, special case of layertype 1
246     newrhoc = rhocobs + porosity*icedensity*cice
247     k = porosity*kice + kobs ! option A, end-member case of type 1, option A
248     !k = kw*kice/((1-porosity)*kice+porosity*kw) ! option B, harmonic average
249     newti = sqrt(newrhoc*k)
250
251  case (5)  ! custom values
252     ! values from Mellon et al. (2004) for ice-cemented soil
253     newrhoc = 2018.*1040.
254     k = 2.5
255     newti = sqrt(newrhoc*k)
256
257  case default
258     error stop 'invalid layer type'
259     
260  end select
261 
262end subroutine soilthprop
263
264
265!=======================================================================
266     
267      real*8 function frostpoint(p)
268!     inverse of psv
269!     input is partial pressure [Pascal]
270!     output is temperature [Kelvin]
271      implicit none
272      real*8 p
273     
274!-----inverse of parametrization 1
275!      real*8 DHmelt,DHvap,DHsub,R,pt,Tt
276!      parameter (DHmelt=6008.,DHvap=45050.)
277!      parameter (DHsub=DHmelt+DHvap)
278!      parameter (R=8.314,pt=6.11e2,Tt=273.16)
279!      frostpoint = 1./(1./Tt-R/DHsub*log(p/pt))
280     
281!-----inverse of parametrization 2
282!     inverse of eq. (2) in Murphy & Koop (2005)
283      real*8 A,B
284      parameter (A=-6143.7, B=28.9074)
285      frostpoint = A / (log(p) - B)
286
287!-----approximate inverse of parametrization 3
288!     eq. (8) in Murphy & Koop (2005)
289!      frostpoint = (1.814625*log(p) + 6190.134)/(29.120 - log(p))
290     
291      end function frostpoint
292
293END MODULE dyn_ss_ice_m_mod
Note: See TracBrowser for help on using the repository browser.