source: trunk/LMDZ.MARS/libf/phymars/phyredem.F90 @ 3026

Last change on this file since 3026 was 2999, checked in by llange, 18 months ago

Mars PCM
Include perenial_co2ice (equivalent of watercap) to distinguich between CO2 frost and perenial CO2 ice for paleoclimate studies.
When no frost is present and we dig into perenial ice, the surface albedo is changed. The albedo for seasonal ice is set to 0.65, and the perenial ice albedo can be fixed in the callphys.def. I recommand values between 0.8 and 0.9.
To use this, paleoclimate must be set to True and TESalbedo to false in the callphys.def. Else, it runs as usual with TES albedo
LL

File size: 12.2 KB
Line 
1module phyredem
2
3implicit none
4
5contains
6
7subroutine physdem0(filename,lonfi,latfi,nsoil,ngrid,nlay,nq, &
8                         phystep,day_ini,time,airefi, &
9                         alb,ith,def_slope, &
10                         subslope_dist)
11! create physics restart file and write time-independent variables
12  use comsoil_h, only: inertiedat, volcapa, mlayer
13  use geometry_mod, only: cell_area
14  use surfdat_h, only: zmea, zstd, zsig, zgam, zthe, &
15                       z0_default, albedice, emisice, emissiv, &
16                       iceradius, dtemisice, phisfi, z0,   &
17                       hmons,summit,base,watercaptag
18  use dimradmars_mod, only: tauvis
19  use iostart, only : open_restartphy, close_restartphy, &
20                      put_var, put_field, length
21  use mod_grid_phy_lmdz, only : klon_glo
22  use planete_h, only: aphelie, emin_turb, lmixmin, obliquit, &
23                       peri_day, periheli, year_day
24  use comcstfi_h, only: g, mugaz, omeg, rad, rcp
25  use time_phylmdz_mod, only: daysec
26  use comslope_mod, ONLY: nslope
27  USE paleoclimate_mod, ONLY: paleoclimate, lag_h2o_ice, lag_co2_ice 
28  implicit none
29 
30  character(len=*), intent(in) :: filename
31  real,intent(in) :: lonfi(ngrid)
32  real,intent(in) :: latfi(ngrid)
33  integer,intent(in) :: nsoil
34  integer,intent(in) :: ngrid
35  integer,intent(in) :: nlay
36  integer,intent(in) :: nq
37  real,intent(in) :: phystep
38  real,intent(in) :: day_ini
39  real,intent(in) :: time
40  real,intent(in) :: airefi(ngrid)
41  real,intent(in) :: alb(ngrid)
42  real,intent(in) :: ith(ngrid,nsoil) ! thermal inertia for present day
43  real, intent(in) :: def_slope(nslope+1) !boundaries for bining of the slopes
44  real, intent(in) :: subslope_dist(ngrid,nslope) !undermesh statistics
45
46  real :: tab_cntrl(length) ! nb "length=100" defined in iostart module
47  integer :: ig
48  real :: watercaptag_tmp(ngrid)
49
50  ! Create physics start file
51  call open_restartphy(filename)
52 
53  ! Build tab_cntrl(:) array
54  tab_cntrl(:)=0.0
55  !cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
56  ! Fill control array tab_cntrl(:) with parameters for this run
57  !cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
58  ! Informations on the physics grid
59  tab_cntrl(1) = float(klon_glo)  ! total number of nodes on physics grid
60  tab_cntrl(2) = float(nlay) ! number of atmospheric layers
61  tab_cntrl(3) = day_ini + int(time)      ! initial day
62  tab_cntrl(4) = time -int(time)          ! initial time of day
63
64  ! Informations about Mars, used by dynamics and physics
65  tab_cntrl(5) = rad      ! radius of Mars (m) ~3397200
66  tab_cntrl(6) = omeg     ! rotation rate (rad.s-1)
67  tab_cntrl(7) = g        ! gravity (m.s-2) ~3.72
68  tab_cntrl(8) = mugaz    ! Molar mass of the atmosphere (g.mol-1) ~43.49
69  tab_cntrl(9) = rcp      !  = r/cp  ~0.256793 (=kappa dans dynamique)
70  tab_cntrl(10) = daysec  ! length of a sol (s)  ~88775
71
72  tab_cntrl(11) = phystep  ! time step in the physics
73  tab_cntrl(12) = 0.
74  tab_cntrl(13) = 0.
75
76  ! Informations about Mars, only for physics
77  tab_cntrl(14) = year_day  ! length of year (sols) ~668.6
78  tab_cntrl(15) = periheli  ! min. Sun-Mars distance (Mkm) ~206.66
79  tab_cntrl(16) = aphelie   ! max. SUn-Mars distance (Mkm) ~249.22
80  tab_cntrl(17) = peri_day  ! date of perihelion (sols since N. spring)
81  tab_cntrl(18) = obliquit  ! Obliquity of the planet (deg) ~23.98
82
83  ! Boundary layer and turbulence
84  tab_cntrl(19) = z0_default   ! default surface roughness (m) ~0.01
85  tab_cntrl(20) = lmixmin   ! mixing length ~100
86  tab_cntrl(21) = emin_turb ! minimal energy ~1.e-8
87
88  ! Optical properties of polar caps and ground emissivity
89  tab_cntrl(22) = albedice(1)  ! Albedo of northern cap ~0.5
90  tab_cntrl(23) = albedice(2)  ! Albedo of southern cap ~0.5
91  tab_cntrl(24) = emisice(1)   ! Emissivity of northern cap ~0.95
92  tab_cntrl(25) = emisice(2)   ! Emissivity of southern cap ~0.95
93  tab_cntrl(26) = emissiv      ! Emissivity of martian soil ~.95
94  tab_cntrl(31) = iceradius(1) ! mean scat radius of CO2 snow (north)
95  tab_cntrl(32) = iceradius(2) ! mean scat radius of CO2 snow (south)
96  tab_cntrl(33) = dtemisice(1) ! time scale for snow metamorphism (north)
97  tab_cntrl(34) = dtemisice(2) ! time scale for snow metamorphism (south)
98
99  ! dust aerosol properties
100  tab_cntrl(27) = tauvis      ! mean visible optical depth
101
102  ! Soil properties:
103  tab_cntrl(35) = volcapa ! soil volumetric heat capacity
104
105  ! Write the controle array
106  call put_var("controle","Control parameters",tab_cntrl)
107 
108  ! Write the mid-layer depths
109  call put_var("soildepth","Soil mid-layer depth",mlayer)
110 
111  ! Write longitudes
112  call put_field("longitude","Longitudes of physics grid",lonfi)
113 
114  ! Write latitudes
115  call put_field("latitude","Latitudes of physics grid",latfi)
116 
117  ! Write mesh areas
118  call put_field("area","Mesh area",cell_area)
119 
120  ! Write surface geopotential
121  call put_field("phisfi","Geopotential at the surface",phisfi)
122 
123  ! Write surface albedo
124  call put_field("albedodat","Albedo of bare ground",alb)
125 
126  ! Subgrid topogaphy variables
127  call put_field("ZMEA","Relief: mean relief",zmea)
128  call put_field("ZSTD","Relief: standard deviation",zstd)
129  call put_field("ZSIG","Relief: sigma parameter",zsig)
130  call put_field("ZGAM","Relief: gamma parameter",zgam)
131  call put_field("ZTHE","Relief: theta parameter",zthe)
132  call put_field("hmons","Relief: hmons parameter (summit - base)",hmons)
133  call put_field("summit","Relief: altitude from the aeroid of the summit of the highest subgrid topography",summit)
134  call put_field("base","Relief: altitude from the aeroid of the base of the highest subgrid topography",base)
135     
136  ! Soil thermal inertia
137  call put_field("inertiedat","Soil thermal inertia - present day",ith)
138 
139  ! Surface roughness length
140  call put_field("z0","Surface roughness length",z0)
141 
142  do ig=1,ngrid
143    if(watercaptag(ig)) then
144       watercaptag_tmp(ig)=1
145    else
146       watercaptag_tmp(ig)=0
147    endif
148  enddo
149
150  call put_field("watercaptag","Infinite water reservoir",watercaptag_tmp)
151
152  ! Sub grid scale slope parametrization
153  call put_var("def_slope","slope criterium stages",def_slope)
154  call put_field("subslope_dist","under mesh slope distribution",subslope_dist)
155
156  ! Paleoclimate outputs
157  if(paleoclimate) then
158    call put_field("lag_h2o_ice","Depth of the H2O lags",lag_h2o_ice)
159    call put_field("lag_co2_ice","Depth of the CO2 lags",lag_co2_ice)
160  endif
161  ! Close file
162  call close_restartphy
163 
164end subroutine physdem0
165
166subroutine physdem1(filename,nsoil,ngrid,nlay,nq, &
167                    phystep,time,tsurf,tsoil,inertiesoil, &
168                    albedo,emis,q2,qsurf,&
169                    tauscaling,totcloudfrac,wstar, &
170                    watercap,perenial_co2ice)
171  ! write time-dependent variable to restart file
172  use iostart, only : open_restartphy, close_restartphy, &
173                      put_var, put_field
174  use tracer_mod, only: noms ! tracer names
175  use nonoro_gwd_ran_mod, only: du_nonoro_gwd, dv_nonoro_gwd
176  use compute_dtau_mod, only: dtau
177  use dust_rad_adjust_mod, only: dust_rad_adjust_prev,dust_rad_adjust_next
178  use dust_param_mod, only: dustscaling_mode
179  use comsoil_h,only: flux_geo
180  use comslope_mod, only: nslope
181  use paleoclimate_mod, only: paleoclimate
182  implicit none
183 
184  include "callkeys.h"
185 
186  character(len=*),intent(in) :: filename
187  integer,intent(in) :: nsoil
188  integer,intent(in) :: ngrid
189  integer,intent(in) :: nlay
190  integer,intent(in) :: nq
191  real,intent(in) :: phystep
192  real,intent(in) :: time
193  real,intent(in) :: tsurf(ngrid,nslope)
194  real,intent(in) :: tsoil(ngrid,nsoil,nslope)
195  real,intent(in) :: inertiesoil(ngrid,nsoil,nslope)
196  real,intent(in) :: albedo(ngrid,2,nslope)
197  real,intent(in) :: emis(ngrid,nslope)
198  real,intent(in) :: q2(ngrid,nlay+1)
199  real,intent(in) :: qsurf(ngrid,nq,nslope)
200  real,intent(in) :: tauscaling(ngrid)
201  real,intent(in) :: totcloudfrac(ngrid)
202  real,intent(in) :: wstar(ngrid)
203  real,intent(in) :: watercap(ngrid,nslope)
204  real,intent(in) :: perenial_co2ice(ngrid,nslope)
205
206  integer :: iq
207  character(len=30) :: txt ! to store some text
208  ! indexes of water vapour & water ice tracers (if any):
209  integer :: i_h2o_vap=0
210  integer :: i_h2o_ice=0
211  integer :: i_hdo_vap=0
212  integer :: i_hdo_ice=0
213
214 
215  ! Open file
216  call open_restartphy(filename)
217 
218  ! First variable to write must be "Time", in order to correctly
219  ! set time counter in file
220  call put_var("Time","Temps de simulation",time)
221
222  ! Water ice layer
223  call put_field("watercap","H2O ice cover",watercap,time)
224 
225 ! Perenial CO2 ice layer
226  if(paleoclimate) then
227    call put_field("perenial_co2ice","CO2 ice cover",perenial_co2ice,time)
228  endif
229  ! Surface temperature
230  call put_field("tsurf","Surface temperature",tsurf,time)
231
232  ! Soil temperature
233  call put_field("inertiesoil","Soil thermal inertia",inertiesoil,time)
234 
235  ! Soil temperature
236  call put_field("tsoil","Soil temperature",tsoil,time)
237 
238  ! Albedo of the surface
239  call put_field("albedo","Surface albedo",albedo(:,1,:),time)
240 
241  ! Emissivity of the surface
242  call put_field("emis","Surface emissivity",emis,time)
243 
244  ! Planetary Boundary Layer
245  call put_field("q2","pbl wind variance",q2,time)
246
247  ! Sub-grid cloud fraction
248  call put_field("totcloudfrac","Total cloud fraction",totcloudfrac,time)
249 
250  ! Dust conversion factor
251  ! Only to be read by newstart to convert to actual dust values
252  ! Or by any user who wants to reconstruct dust, opacity from the start files.
253  call put_field("tauscaling","dust conversion factor",tauscaling,time)
254
255  ! Radiative scaling coefficients
256  if (dustscaling_mode==2) then
257    call put_field("dust_rad_adjust_prev", &
258                   "radiative dust adjustement factor prev. sol", &
259                   dust_rad_adjust_prev,time)
260    call put_field("dust_rad_adjust_next", &
261                   "radiative dust adjustement factor next sol", &
262                   dust_rad_adjust_next,time)
263  endif
264 
265  if (dustinjection.gt.0) then
266    call put_field("dtau","dust opacity difference between GCM and scenario",&
267                   dtau,time)
268  endif
269
270  if (calltherm) then
271    call put_field("wstar","Max vertical velocity in thermals",wstar,time)
272  endif
273
274  ! Tracers on the surface
275  ! preliminary stuff: look for water vapour & water ice tracers (if any)
276  do iq=1,nq
277    if (noms(iq).eq."h2o_vap") then
278      i_h2o_vap=iq
279    endif
280    if (noms(iq).eq."h2o_ice") then
281      i_h2o_ice=iq
282    endif
283
284  ! look for HDO vapour & HDO ice tracers (if any)
285    if (noms(iq).eq."hdo_vap") then
286      i_hdo_vap=iq
287    endif
288    if (noms(iq).eq."hdo_ice") then
289      i_hdo_ice=iq
290    endif
291  enddo
292
293 
294  if (nq.gt.0) then
295    do iq=1,nq
296      txt=noms(iq)
297      ! Exception: there is no water vapour surface tracer
298      if (txt.eq."h2o_vap") then
299        write(*,*)"physdem1: skipping water vapour tracer"
300        if (i_h2o_ice.eq.i_h2o_vap) then
301          ! then there is no "water ice" tracer; but still
302          ! there is some water ice on the surface
303          write(*,*)"          writing water ice instead"
304          txt="h2o_ice"
305        else
306          ! there is a "water ice" tracer which has been / will be
307          ! delt with in due time
308          cycle
309        endif ! of if (igcm_h2o_ice.eq.igcm_h2o_vap)
310      endif ! of if (txt.eq."h2o_vap")
311
312      if (txt.eq."hdo_vap") then
313        write(*,*)"physdem1: skipping HDO vapour tracer"
314        if (i_hdo_ice.eq.i_hdo_vap) then
315          ! then there is no "water ice" tracer; but still
316          ! there is some water ice on the surface
317          write(*,*)"          writing HDO ice instead"
318          txt="hdo_ice"
319        else
320          ! there is a "water ice" tracer which has been / will be
321          ! delt with in due time
322          cycle
323        endif ! of if (igcm_hdo_ice.eq.igcm_hdo_vap)
324      endif ! of if (txt.eq."hdo_vap")
325
326      ! co2_ice has been added to qsurf(:,igcm_co2) in co2condens4micro
327      if (txt.eq."co2_ice") then
328        write(*,*)"physdem1: skipping co2_ice tracer"
329        cycle
330      end if     
331
332      call put_field(trim(txt),"tracer on surface",qsurf(:,iq,:),time)
333    enddo
334  endif
335
336  ! Non-orographic gavity waves
337  if (calllott_nonoro) then
338     call put_field("du_nonoro_gwd","Zonal wind tendency due to GW",du_nonoro_gwd,time)
339     call put_field("dv_nonoro_gwd","Meridional wind tendency due to GW",dv_nonoro_gwd,time)
340  endif
341
342
343  ! Geothermal Flux
344     call put_field('flux_geo','Geothermal flux',flux_geo,time)
345  ! Close file
346  call close_restartphy
347 
348end subroutine physdem1
349
350end module phyredem
Note: See TracBrowser for help on using the repository browser.