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

Last change on this file since 2932 was 2913, checked in by romain.vande, 22 months ago

Mars PCM:
Adapt start2archive.F to the subslope parametrisation.
Small correction for some dimensions of variables.
RV

File size: 11.5 KB
RevLine 
[1130]1module phyredem
2
3implicit none
4
5contains
6
7subroutine physdem0(filename,lonfi,latfi,nsoil,ngrid,nlay,nq, &
8                         phystep,day_ini,time,airefi, &
[2896]9                         alb,ith,def_slope, &
10                         subslope_dist)
[1130]11! create physics restart file and write time-independent variables
12  use comsoil_h, only: inertiedat, volcapa, mlayer
[1543]13  use geometry_mod, only: cell_area
[1221]14  use surfdat_h, only: zmea, zstd, zsig, zgam, zthe, &
[1130]15                       z0_default, albedice, emisice, emissiv, &
[1974]16                       iceradius, dtemisice, phisfi, z0,   &
[2884]17                       hmons,summit,base,watercaptag
[1246]18  use dimradmars_mod, only: tauvis
[1130]19  use iostart, only : open_restartphy, close_restartphy, &
20                      put_var, put_field, length
21  use mod_grid_phy_lmdz, only : klon_glo
[1524]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
[2896]26  use comslope_mod, ONLY: nslope
[1130]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)
[2896]42  real, intent(in) :: def_slope(nslope+1) !boundaries for bining of the slopes
43  real, intent(in) :: subslope_dist(ngrid,nslope) !undermesh statistics
[1974]44
[1130]45  real :: tab_cntrl(length) ! nb "length=100" defined in iostart module
[2884]46  integer :: ig
47  real :: watercaptag_tmp(ngrid)
[1130]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
[1541]117  call put_field("area","Mesh area",cell_area)
[1130]118 
119  ! Write surface geopotential
120  call put_field("phisfi","Geopotential at the surface",phisfi)
121 
122  ! Write surface albedo
[1221]123  call put_field("albedodat","Albedo of bare ground",alb)
[1130]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)
[2079]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)
[1974]134     
[1130]135  ! Soil thermal inertia
[1210]136  call put_field("inertiedat","Soil thermal inertia",ith)
[1130]137 
138  ! Surface roughness length
139  call put_field("z0","Surface roughness length",z0)
140 
[2884]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
[2896]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
[1130]155  ! Close file
156  call close_restartphy
157 
158end subroutine physdem0
159
160subroutine physdem1(filename,nsoil,ngrid,nlay,nq, &
[2826]161                    phystep,time,tsurf,tsoil,albedo,emis,q2,qsurf,&
[1944]162                    tauscaling,totcloudfrac,wstar, &
[2562]163                    watercap)
[1130]164  ! write time-dependent variable to restart file
165  use iostart, only : open_restartphy, close_restartphy, &
166                      put_var, put_field
[1621]167  use tracer_mod, only: noms ! tracer names
[2220]168  use nonoro_gwd_ran_mod, only: du_nonoro_gwd, dv_nonoro_gwd
[2265]169  use compute_dtau_mod, only: dtau
[2417]170  use dust_rad_adjust_mod, only: dust_rad_adjust_prev,dust_rad_adjust_next
171  use dust_param_mod, only: dustscaling_mode
[2887]172  use comsoil_h,only: flux_geo
[2913]173  use comslope_mod, ONLY: nslope
[1130]174  implicit none
[1922]175 
176  include "callkeys.h"
177 
[1130]178  character(len=*),intent(in) :: filename
179  integer,intent(in) :: nsoil
180  integer,intent(in) :: ngrid
181  integer,intent(in) :: nlay
182  integer,intent(in) :: nq
183  real,intent(in) :: phystep
184  real,intent(in) :: time
[2913]185  real,intent(in) :: tsurf(ngrid,nslope)
186  real,intent(in) :: tsoil(ngrid,nsoil,nslope)
187  real,intent(in) :: albedo(ngrid,2,nslope)
188  real,intent(in) :: emis(ngrid,nslope)
[1130]189  real,intent(in) :: q2(ngrid,nlay+1)
[2913]190  real,intent(in) :: qsurf(ngrid,nq,nslope)
[1208]191  real,intent(in) :: tauscaling(ngrid)
[1922]192  real,intent(in) :: totcloudfrac(ngrid)
[1944]193  real,intent(in) :: wstar(ngrid)
[2913]194  real,intent(in) :: watercap(ngrid,nslope)
[1130]195 
196  integer :: iq
197  character(len=30) :: txt ! to store some text
198  ! indexes of water vapour & water ice tracers (if any):
199  integer :: i_h2o_vap=0
200  integer :: i_h2o_ice=0
[2312]201  integer :: i_hdo_vap=0
202  integer :: i_hdo_ice=0
[1130]203
204 
205  ! Open file
206  call open_restartphy(filename)
207 
208  ! First variable to write must be "Time", in order to correctly
209  ! set time counter in file
210  call put_var("Time","Temps de simulation",time)
[2260]211
212  ! Water ice layer
213  call put_field("watercap","H2O ice cover",watercap,time)
[1130]214 
215  ! Surface temperature
216  call put_field("tsurf","Surface temperature",tsurf,time)
217 
218  ! Soil temperature
219  call put_field("tsoil","Soil temperature",tsoil,time)
220 
[1944]221  ! Albedo of the surface
[2913]222  call put_field("albedo","Surface albedo",albedo(:,1,:),time)
[1944]223 
[1130]224  ! Emissivity of the surface
225  call put_field("emis","Surface emissivity",emis,time)
226 
227  ! Planetary Boundary Layer
228  call put_field("q2","pbl wind variance",q2,time)
[1711]229
230  ! Sub-grid cloud fraction
231  call put_field("totcloudfrac","Total cloud fraction",totcloudfrac,time)
232 
[1208]233  ! Dust conversion factor
234  ! Only to be read by newstart to convert to actual dust values
235  ! Or by any user who wants to reconstruct dust, opacity from the start files.
236  call put_field("tauscaling","dust conversion factor",tauscaling,time)
237
[2417]238  ! Radiative scaling coefficients
239  if (dustscaling_mode==2) then
240    call put_field("dust_rad_adjust_prev", &
241                   "radiative dust adjustement factor prev. sol", &
242                   dust_rad_adjust_prev,time)
243    call put_field("dust_rad_adjust_next", &
244                   "radiative dust adjustement factor next sol", &
245                   dust_rad_adjust_next,time)
246  endif
247 
[2265]248  if (dustinjection.gt.0) then
249    call put_field("dtau","dust opacity difference between GCM and scenario",&
250                   dtau,time)
251  endif
252
[1944]253  if (calltherm) then
254    call put_field("wstar","Max vertical velocity in thermals",wstar,time)
255  endif
256
[1130]257  ! Tracers on the surface
258  ! preliminary stuff: look for water vapour & water ice tracers (if any)
259  do iq=1,nq
[1621]260    if (noms(iq).eq."h2o_vap") then
[1130]261      i_h2o_vap=iq
262    endif
[1621]263    if (noms(iq).eq."h2o_ice") then
[1130]264      i_h2o_ice=iq
265    endif
[2312]266
267  ! look for HDO vapour & HDO ice tracers (if any)
268    if (noms(iq).eq."hdo_vap") then
269      i_hdo_vap=iq
270    endif
271    if (noms(iq).eq."hdo_ice") then
272      i_hdo_ice=iq
273    endif
[1130]274  enddo
[2312]275
[1130]276 
277  if (nq.gt.0) then
278    do iq=1,nq
[1621]279      txt=noms(iq)
[1130]280      ! Exception: there is no water vapour surface tracer
281      if (txt.eq."h2o_vap") then
282        write(*,*)"physdem1: skipping water vapour tracer"
283        if (i_h2o_ice.eq.i_h2o_vap) then
284          ! then there is no "water ice" tracer; but still
285          ! there is some water ice on the surface
286          write(*,*)"          writing water ice instead"
287          txt="h2o_ice"
288        else
289          ! there is a "water ice" tracer which has been / will be
290          ! delt with in due time
291          cycle
292        endif ! of if (igcm_h2o_ice.eq.igcm_h2o_vap)
293      endif ! of if (txt.eq."h2o_vap")
[2312]294
295      if (txt.eq."hdo_vap") then
296        write(*,*)"physdem1: skipping HDO vapour tracer"
297        if (i_hdo_ice.eq.i_hdo_vap) then
298          ! then there is no "water ice" tracer; but still
299          ! there is some water ice on the surface
300          write(*,*)"          writing HDO ice instead"
301          txt="hdo_ice"
302        else
303          ! there is a "water ice" tracer which has been / will be
304          ! delt with in due time
305          cycle
306        endif ! of if (igcm_hdo_ice.eq.igcm_hdo_vap)
307      endif ! of if (txt.eq."hdo_vap")
308
[2826]309      ! co2_ice has been added to qsurf(:,igcm_co2) in co2condens4micro
[2494]310      if (txt.eq."co2_ice") then
311        write(*,*)"physdem1: skipping co2_ice tracer"
312        cycle
313      end if     
314
[2913]315      call put_field(trim(txt),"tracer on surface",qsurf(:,iq,:),time)
[1130]316    enddo
317  endif
[2562]318
[2220]319  ! Non-orographic gavity waves
320  if (calllott_nonoro) then
321     call put_field("du_nonoro_gwd","Zonal wind tendency due to GW",du_nonoro_gwd,time)
322     call put_field("dv_nonoro_gwd","Meridional wind tendency due to GW",dv_nonoro_gwd,time)
323  endif
[2494]324
[2887]325
326  ! Geothermal Flux
327     call put_field('flux_geo','Geothermal flux',flux_geo,time)
[1130]328  ! Close file
329  call close_restartphy
330 
331end subroutine physdem1
332
333end module phyredem
Note: See TracBrowser for help on using the repository browser.