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

Last change on this file since 3495 was 3493, checked in by jbclement, 3 weeks ago

PEM:

  • Renaming of Norbert Schorghofer's subroutines with the prefix 'NS_';
  • Making the extension of all NS's subroutines as '.F90';
  • Deletion of the wrapper subroutine;
  • Making the initialization, variables management and arguments of the main subroutine for the dynamic computation of ice table to be more suitable.

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