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

Last change on this file since 3493 was 3338, checked in by jbclement, 8 months ago

Mars PCM:
Addition of the paleoclimate variable 'd_coef' in the writing of the "restartfi.nc" file.
JBC

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