source: trunk/LMDZ.COMMON/libf/evolution/pemredem.F90 @ 2850

Last change on this file since 2850 was 2842, checked in by romain.vande, 3 years ago

Mars PEM:
PEM is now adapted to run with XIOS diurnal averages (when they will work properly)
RV

File size: 4.4 KB
Line 
1module pemredem
2
3implicit none
4
5contains
6
7subroutine pemdem0(filename,lonfi,latfi,cell_area,nsoil_PEM,ngrid, &
8                         day_ini,time,nslope,def_slope,subslope_dist)
9! create physics restart file and write time-independent variables
10  use comsoil_h_PEM, only: mlayer_PEM,fluxgeo,inertiedat_PEM
11  use iostart_PEM, only : open_restartphy, close_restartphy, &
12                      put_var, put_field, length
13  use mod_grid_phy_lmdz, only : klon_glo
14#ifndef CPP_STD
15  use planete_h, only: aphelie, emin_turb, lmixmin, obliquit, &
16                       peri_day, periheli, year_day
17  use comcstfi_h, only: g, mugaz, omeg, rad, rcp
18#else
19  use planete_mod, only: apoastr, emin_turb, lmixmin, obliquit, &
20                       peri_day, periastr, year_day
21  USE comconst_mod, ONLY: g, omeg, rad
22  use comcstfi_mod, only: mugaz, rcp
23#endif
24  use time_phylmdz_mod, only: daysec
25  implicit none
26 
27  character(len=*), intent(in) :: filename
28  real,intent(in) :: lonfi(ngrid)
29  real,intent(in) :: latfi(ngrid)
30  integer,intent(in) :: nsoil_PEM
31  integer,intent(in) :: ngrid
32  real,intent(in) :: day_ini
33  real,intent(in) :: time
34  real, intent(in) :: cell_area(ngrid) !boundaries for bining of the slopes
35  integer,intent(in) :: nslope
36  real, intent(in) :: def_slope(nslope+1) !boundaries for bining of the slopes
37  real, intent(in) :: subslope_dist(ngrid,nslope) !undermesh statistics
38 
39  ! Create physics start file
40  call open_restartphy(filename)
41 
42  ! Write the mid-layer depths
43  call put_var("soildepth","Soil mid-layer depth",mlayer_PEM)
44 
45  ! Write longitudes
46  call put_field("longitude","Longitudes of physics grid",lonfi)
47 
48  ! Write latitudes
49  call put_field("latitude","Latitudes of physics grid",latfi)
50 
51  ! Write mesh areas
52  call put_field("area","Mesh area",cell_area)
53     
54  ! Multidimensionnal variables (nogcm undermesh slope statistics)
55  call put_var("def_slope","slope criterium stages",def_slope)
56  call put_field("subslope_dist","under mesh slope distribution",subslope_dist)
57
58  call put_var("FluxGeo","Geothermal Flux (W/m^2)",fluxgeo)
59
60  ! Close file
61  call close_restartphy
62 
63end subroutine pemdem0
64
65subroutine pemdem1(filename,year_iter,nsoil_PEM,ngrid,nslope, &
66                    tsoil_slope_PEM,inertiesoil_slope_PEM,ice_table,m_co2_regolith)
67  ! write time-dependent variable to restart file
68  use iostart_PEM, only : open_restartphy, close_restartphy, &
69                      put_var, put_field
70
71  use comsoil_h_PEM, only: mlayer_PEM,fluxgeo,inertiedat_PEM
72  USE temps_mod_evol, ONLY: year_PEM
73  use soil_evolution_mod, only: soil_pem
74               
75  implicit none
76 
77#ifndef CPP_STD
78  include "callkeys.h"
79#endif
80 
81  character(len=*),intent(in) :: filename
82  integer,intent(in) :: nsoil_PEM
83  integer,intent(in) :: ngrid
84  real,intent(IN) :: tsoil_slope_PEM(ngrid,nsoil_PEM,nslope) !under mesh bining according to slope
85  real,intent(IN) :: inertiesoil_slope_PEM(ngrid,nsoil_PEM,nslope) !under mesh bining according to slope
86  real,intent(IN) :: ice_table(ngrid,nslope) !under mesh bining according to slope
87  integer,intent(IN) :: nslope
88  real,intent(in) :: m_co2_regolith(ngrid,nsoil_PEM,nslope)
89  integer :: islope
90  CHARACTER*2 :: num 
91  integer, intent(IN) :: year_iter
92
93  real :: year_tot
94 
95  integer :: iq
96  character(len=30) :: txt ! to store some text
97  ! indexes of water vapour & water ice tracers (if any):
98
99  year_tot=real(year_iter+year_PEM)
100  print *, "Year of simulation=", year_tot 
101
102  ! Open file
103  call open_restartphy(filename)
104 
105  ! First variable to write must be "Time", in order to correctly
106  ! set time counter in file
107  call put_var("Time","Year of simulation",year_tot)
108
109    if(soil_pem) then
110 
111  ! Multidimensionnal variables (undermesh slope statistics)
112DO islope=1,nslope
113  write(num,fmt='(i2.2)') islope
114  call put_field("tsoil_PEM_slope"//num,"Soil temperature by slope type", &
115                 tsoil_slope_PEM(:,:,islope),year_tot)
116  call put_field("TI_PEM_slope"//num,"Soil Thermal Inertia by slope type", &
117                 inertiesoil_slope_PEM(:,:,islope),year_tot)
118
119  call put_field("mco2_reg_ads_slope"//num, "Mass of co2 adsorbded in the regolith", &
120                    m_co2_regolith(:,:,islope), year_tot)
121
122ENDDO
123
124  call put_field("ice_table","Depth of ice table", &
125                 ice_table,year_tot)
126
127  call put_field("inertiedat_PEM","Thermal inertie of PEM ", &
128                 inertiedat_PEM,year_tot)
129
130   endif !soil_pem
131
132  ! Close file
133  call close_restartphy
134 
135end subroutine pemdem1
136
137end module pemredem
Note: See TracBrowser for help on using the repository browser.