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

Last change on this file since 1130 was 1130, checked in by emillour, 11 years ago

Mars GCM:
Series of changes to enable running in parallel (using LMDZ.COMMON dynamics);
Current LMDZ.MARS can still notheless be compiled and run in serial mode
"as previously".
Summary of main changes:

  • Main programs (newstart, start2archive, xvik) that used to be in dyn3d have been moved to phymars.
  • dyn3d/control.h is now module control_mod.F90
  • rearanged input/outputs routines everywhere to handle serial/MPI cases. physdem.F => phyredem.F90 , phyetat0.F => phyetat0.F90 ; all read/write routines for startfi files are gathered in module iostart.F90
  • added parallelism related routines init_phys_lmdz.F90, comgeomphy.F90, dimphy.F90, iniphysiq.F90, mod_grid_phy_lmdz.F90, mod_phys_lmdz_mpi_data.F90, mod_phys_lmdz_mpi_transfert.F90, mod_phys_lmdz_omp_data.F90, mod_phys_lmdz_omp_transfert.F90, mod_phys_lmdz_para.F90, mod_phys_lmdz_transfert_para.F90 in phymars and mod_const_mpi.F90 in dyn3d (for compliance with parallel case)
  • created generic routines 'planetwide_maxval' and 'planetwide_minval', in module "planetwide_mod", that enable obtaining the min and max of a field over the whole planet.

EM

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