1 | module 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 | |
---|
11 | implicit none |
---|
12 | |
---|
13 | contains |
---|
14 | |
---|
15 | subroutine 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 comcstfi_mod, only: g, mugaz, omeg, rad, rcp |
---|
27 | #endif |
---|
28 | |
---|
29 | implicit none |
---|
30 | |
---|
31 | character(len=*), intent(in) :: filename |
---|
32 | real, intent(in) :: lonfi(ngrid) |
---|
33 | real, intent(in) :: latfi(ngrid) |
---|
34 | integer, intent(in) :: nsoil_PEM |
---|
35 | integer, intent(in) :: ngrid |
---|
36 | real, intent(in) :: day_ini |
---|
37 | real, intent(in) :: time |
---|
38 | real, intent(in) :: cell_area(ngrid) !boundaries for bining of the slopes |
---|
39 | integer, intent(in) :: nslope |
---|
40 | real, intent(in) :: def_slope(nslope+1) !boundaries for bining of the slopes |
---|
41 | real, intent(in) :: subslope_dist(ngrid,nslope) !undermesh statistics |
---|
42 | |
---|
43 | ! Create physics start file |
---|
44 | call open_restartphy(filename) |
---|
45 | |
---|
46 | ! Write the mid-layer depths |
---|
47 | call put_var("soildepth","Soil mid-layer depth",mlayer_PEM) |
---|
48 | |
---|
49 | ! Write longitudes |
---|
50 | call put_field("longitude","Longitudes of physics grid",lonfi) |
---|
51 | |
---|
52 | ! Write latitudes |
---|
53 | call put_field("latitude","Latitudes of physics grid",latfi) |
---|
54 | |
---|
55 | ! Write mesh areas |
---|
56 | call put_field("area","Mesh area",cell_area) |
---|
57 | |
---|
58 | ! Multidimensionnal variables (nogcm undermesh slope statistics) |
---|
59 | call put_var("def_slope","slope criterium stages",def_slope) |
---|
60 | call put_field("subslope_dist","under mesh slope distribution",subslope_dist) |
---|
61 | |
---|
62 | |
---|
63 | ! Close file |
---|
64 | call close_restartphy |
---|
65 | |
---|
66 | end subroutine pemdem0 |
---|
67 | |
---|
68 | subroutine pemdem1(filename,i_myear,nsoil_PEM,ngrid,nslope, & |
---|
69 | tsoil_slope_PEM,inertiesoil_slope_PEM,ice_table_depth,ice_table_thickness, & |
---|
70 | m_co2_regolith,m_h2o_regolith,water_reservoir) |
---|
71 | |
---|
72 | |
---|
73 | ! write time-dependent variable to restart file |
---|
74 | use iostart_PEM, only: open_restartphy, close_restartphy, put_var, put_field |
---|
75 | use comsoil_h_PEM, only: mlayer_PEM,fluxgeo, inertiedat_PEM,soil_pem |
---|
76 | use time_evol_mod, only: year_bp_ini, convert_years |
---|
77 | |
---|
78 | implicit none |
---|
79 | |
---|
80 | #ifndef CPP_STD |
---|
81 | include "callkeys.h" |
---|
82 | #endif |
---|
83 | |
---|
84 | character(len=*), intent(in) :: filename |
---|
85 | integer, intent(in) :: nsoil_PEM, ngrid, nslope, i_myear |
---|
86 | real, intent(in) :: tsoil_slope_PEM(ngrid,nsoil_PEM,nslope) !under mesh bining according to slope |
---|
87 | real, intent(in) :: inertiesoil_slope_PEM(ngrid,nsoil_PEM,nslope) !under mesh bining according to slope |
---|
88 | real, intent(in) :: ice_table_depth(ngrid,nslope) !under mesh bining according to slope |
---|
89 | real, intent(in) :: ice_table_thickness(ngrid,nslope) !under mesh bining according to slope |
---|
90 | real, intent(in) :: m_co2_regolith(ngrid,nsoil_PEM,nslope) |
---|
91 | real, intent(in) :: m_h2o_regolith(ngrid,nsoil_PEM,nslope) |
---|
92 | real, intent(in) :: water_reservoir(ngrid) |
---|
93 | |
---|
94 | integer :: islope |
---|
95 | character*2 :: num |
---|
96 | integer :: iq |
---|
97 | character(len=30) :: txt ! To store some text |
---|
98 | real :: Year ! Year of the simulation |
---|
99 | |
---|
100 | ! Open file |
---|
101 | call open_restartphy(filename) |
---|
102 | |
---|
103 | ! First variable to write must be "Time", in order to correctly |
---|
104 | ! set time counter in file |
---|
105 | Year = (year_bp_ini + i_myear)*convert_years |
---|
106 | call put_var("Time","Year of simulation",Year) |
---|
107 | call put_field('water_reservoir','water_reservoir',water_reservoir,Year) |
---|
108 | |
---|
109 | if(soil_pem) then |
---|
110 | |
---|
111 | ! Multidimensionnal variables (undermesh slope statistics) |
---|
112 | do 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) |
---|
116 | call put_field("TI_PEM_slope"//num,"Soil Thermal Inertia by slope type", & |
---|
117 | inertiesoil_slope_PEM(:,:,islope),Year) |
---|
118 | |
---|
119 | call put_field("mco2_reg_ads_slope"//num, "Mass of co2 adsorbded in the regolith", & |
---|
120 | m_co2_regolith(:,:,islope),Year) |
---|
121 | call put_field("mh2o_reg_ads_slope"//num, "Mass of h2o adsorbded in the regolith", & |
---|
122 | m_h2o_regolith(:,:,islope),Year) |
---|
123 | enddo |
---|
124 | |
---|
125 | call put_field("ice_table_depth","Depth of ice table",ice_table_depth,Year) |
---|
126 | |
---|
127 | call put_field("ice_table_thickness","Depth of ice table",ice_table_thickness,Year) |
---|
128 | |
---|
129 | call put_field("inertiedat_PEM","Thermal inertie of PEM ",inertiedat_PEM,Year) |
---|
130 | |
---|
131 | endif !soil_pem |
---|
132 | |
---|
133 | ! Close file |
---|
134 | call close_restartphy |
---|
135 | |
---|
136 | end subroutine pemdem1 |
---|
137 | |
---|
138 | end module pemredem |
---|