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

Last change on this file since 2616 was 2562, checked in by cmathe, 3 years ago

GCM MARS: CO2 clouds microphysics improvements

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