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

Last change on this file since 3512 was 3512, checked in by jbclement, 9 days ago

PEM:
Few corrections related to r3498 (time step from integer to real) and r3493 (Norbert Schorghofer's subroutines for dynamic ice table) in order to make the code work properly.
JBC

File size: 5.8 KB
Line 
1
2!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
3!!!
4!!! Purpose:  Retreat and growth of subsurface ice on Mars
5!!!           orbital elements remain constant
6!!!
7!!!
8!!! Author: EV, updated NS MSIM dynamical program for the PEM
9!!!
10!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
11
12SUBROUTINE dyn_ss_ice_m(ssi_depth_in,T1,Tb,nz,thIn,p0,pfrost,porefill_in,porefill,ssi_depth)
13
14!***********************************************************************
15! Retreat and growth of subsurface ice on Mars
16! orbital elements remain constant
17!***********************************************************************
18  use miscparameters, only : pi, d2r, NMAX, marsDay, solsperyear
19  !use allinterfaces
20  implicit none
21  integer, parameter :: NP=1   ! # of sites
22  integer nz, i, k, iloop
23  real(8)  zmax, delta, z(NMAX), icetime, porosity, icefrac
24  real(8), dimension(NP) ::  albedo, thIn, rhoc
25  real(8), dimension(NP) :: pfrost, p0
26  real(8) newti, stretch, newrhoc, ecc, omega, eps, timestep
27  real(8)  ssi_depth_in, ssi_depth, T1
28  real(8), dimension(NP) :: zdepthF, zdepthE, zdepthT, zdepthG
29  real(8), dimension(NMAX,NP) :: porefill, porefill_in
30  real(8), dimension(nz) :: Tb
31  real(8), dimension(NP) :: Tmean1, Tmean3, avrho1
32  real(8) tmax, tlast, avrho1prescribed(NP), l1
33  real(8), external :: smartzfac
34
35  !if (iargc() /= 1) then
36  !   stop 'USAGE: icages ext'
37  !endif
38  !call getarg( 1, ext )
39
40  if (NP>100) stop 'subroutine icelayer_mars cannot handle this many sites'
41
42  ! parameters that never ever change
43  porosity = 0.4d0  ! porosity of till
44  !rhoc(:) = 1500.*800.  ! will be overwritten
45  icefrac = 0.98
46  tmax = 1
47  tlast = 0.
48  avrho1prescribed(:) = pfrost/T1  ! <0 means absent
49  albedo=0.23
50  !avrho1prescribed(:) = 0.16/200.  ! units are Pa/K
51
52  !open(unit=21,file='lats.'//ext,action='read',status='old',iostat=ierr)
53  !if (ierr /= 0) then
54  !   print *,'File lats.'//ext,'not found'
55  !   stop
56  !endif
57  do k=1,NP
58     !read(21,*) latitude(k),albedo(k),thIn(k),htopo(k)
59     ! empirical relation from Mellon & Jakosky
60     rhoc(k) = 800.*(150.+100.*sqrt(34.2+0.714*thIn(k)))
61  enddo
62  !close(21)
63
64  ! set eternal grid
65  zmax = 25.
66  !zfac = smartzfac(nz,zmax,6,0.032d0)
67  !call setgrid(nz,z,zmax,zfac)
68  l1=2.e-4
69  do iloop=0,nz - 1
70    z(iloop + 1) = l1*(1+iloop**2.9*(1-exp(-real(iloop)/20.)))
71  enddo
72 
73
74  !open(unit=30,file='z.'//ext,action='write',status='unknown')
75  !write(30,'(999(f8.5,1x))') z(1:nz)
76  !close(30)
77
78  !ecc = ecc_in;  eps = obl_in*d2r;  omega = Lp_in*d2r   ! today
79  ! total atmospheric pressure
80  !p0(:) = 600.
81  ! presently 520 Pa at zero elevation (Smith & Zuber, 1998)
82 ! do k=1,NP
83 !    p0(k)=520*exp(-htopo(k)/10800.)
84 ! enddo
85  timestep = 1 ! must be integer fraction of 1 ka
86  icetime = -tmax-timestep  ! earth years
87 
88  ! initializations
89  !Tb = -9999.
90  zdepthF(:) = -9999.
91
92  !zdepthT(1:NP) = -9999.   ! reset again below
93 ! zdepthT(1:NP) = 0.
94
95 ! print *,'RUNNING MARS_FAST'
96 ! print *,'Global model parameters:'
97 ! print *,'nz=',nz,' zfac=',zfac,'zmax=',zmax
98 ! print *,'porosity=',porosity
99 ! print *,'starting at time',icetime,'years'
100 ! print *,'time step=',timestep,'years'
101 ! print *,'eps=',eps/d2r,'ecc=',ecc,'omega=',omega/d2r
102 ! print *,'number of sites=',NP
103 ! print *,'Site specific parameters:'
104  do k=1,NP
105     if (NP>1) print *,'  Site ',k
106 !    print *,'  latitude (deg)',latitude(k),' rho*c (J/m^3/K)',rhoc(k),' thIn=',thIn(k)
107 !    print *,'  total pressure=',p0(k),'partial pressure=',pfrost(k)
108     delta = thIn(k)/rhoc(k)*sqrt(marsDay/pi)
109 !    print *,'  skin depths (m)',delta,delta*sqrt(solsperyear)
110     call soilthprop(porosity,1.d0,rhoc(k),thIn(k),1,newrhoc,newti,icefrac)
111     stretch = (newti/thIn(k))*(rhoc(k)/newrhoc)
112     do i=1,nz
113        if (z(i)<delta) cycle
114 !       print *,'  ',i-1,' grid points within diurnal skin depth'
115        exit
116     enddo
117 !    print *,'  ',zmax/(sqrt(solsperyear)*delta),'times seasonal dry skin depth'
118 !    print *,'  ',zmax/(sqrt(solsperyear)*delta*stretch),'times seasonal filled skin depth'
119 !    print *,'  Initial ice depth=',zdepthT(k)
120 !    print *
121  enddo
122 ! call outputmoduleparameters
123 ! print *
124
125  ! open and name all output files
126!  open(unit=34,file='subout.'//ext,action='write',status='unknown')
127!  open(unit=36,file='depthF.'//ext,action='write',status='unknown')
128!  open(unit=37,file='depths.'//ext,action='write',status='unknown')
129
130 ! print *,'Equilibrating initial temperature'
131 ! do i=1,4
132 !    call icelayer_mars(0d0,nz,NP,thIn,rhoc,z,porosity,pfrost,Tb,zdepthF, &
133 !      &  zdepthE,porefill(1:nz,:),Tmean1,Tmean3,zdepthG, &
134 !      &  latitude,albedo,p0,ecc,omega,eps,icefrac,zdepthT,avrho1, &
135 !      &  avrho1prescribed)
136 ! enddo
137
138  !print *,'History begins here'
139 porefill(1:nz,1:NP) =  porefill_in(1:nz,1:NP)
140  zdepthT(1:NP) = ssi_depth_in
141  do
142  !print *,'Zt0=  ',ZdepthT
143     call icelayer_mars(timestep,nz,NP,thIn,rhoc,z,porosity,pfrost,Tb,zdepthF, &
144          & zdepthE,porefill(1:nz,:),Tmean1,Tmean3,zdepthG, &
145          & albedo,p0,icefrac,zdepthT,avrho1, &
146          & avrho1prescribed)
147     icetime = icetime+timestep
148   !  print *,'T_after=  ',Tb(:)
149 !    print *,'z=  ',z(:)
150 !    print *,'Zt=  ',ZdepthT
151     ssi_depth=ZdepthT(1)
152    ! if (abs(mod(icetime/100.,1.d0))<1.e-3) then ! output every 1000 years
153    !    do k=1,NP
154         !write(36,*) icetime,latitude(k),zdepthF(k),porefill(1:nz,k)
155           ! compact output format
156 !          write(36,'(f10.0,2x,f7.3,1x,f11.5,1x)',advance='no') &
157 !               & icetime,latitude(k),zdepthF(k)
158 !          call compactoutput(36,porefill(:,k),nz)
159 !          write(37,501) icetime,latitude(k),zdepthT(k), &
160 !               & Tmean1(k),Tmean3(k),zdepthG(k),avrho1(k)
161  !      enddo
162  !   endif
163!     print *,icetime
164     if (icetime>=tlast) exit
165  enddo
166
167 ! close(34)
168 ! close(36); close(37)
169
170!501 format (f10.0,2x,f7.3,2x,f10.4,2(2x,f6.2),2x,f9.3,2x,g11.4)
171 
172end subroutine dyn_ss_ice_m
173
174
175
Note: See TracBrowser for help on using the repository browser.