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

Last change on this file since 2156 was 2079, checked in by mvals, 6 years ago

Mars GCM:

  • Update of rocketduststorm_mod.F90 :

We want to separate both parametrizations related to the formation of the detached dust layers. Therefore, rocketduststorm_mod.F90 now only comprises the rocket dust storm scheme, whereas it contained
also before the calculation of the vertical velocity induced by the presence of the sub-grid scale topography. This latter part is under development and will be integrated as a fully independant
parametrization: the aim is to simulate the entrainment of dust by slope winds, from the boundary layer up to the top of the sub-grid scale topography.

  • Addition of initial surface parameters "summit" and "base" to prepare the previously described slope wind parametrization

MV

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