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

Last change on this file since 2947 was 2942, checked in by llange, 21 months ago

Mars PCM
Add the possibility to use a different thermal inertia (field
'inertiesoil') than inertiedat in the PCM (for paleoclimate studies). By defaut, if there is not
inertiesoil, inertiedat is used. Soil_tifeedback still work with
inertiedat
Newstart adapted, start2archive will be modified by Romain
Vandemeulebrouck.
LL

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