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

Last change on this file since 3040 was 3039, checked in by jbclement, 2 years ago

Mars PEM:
Big cleaning of main program pem.F90 (indentation, declarations, comments, simplification of conditions/loops, etc).
JBC

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