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

Last change on this file since 3552 was 3537, checked in by jbclement, 2 weeks ago

PEM:
Small fixes to initialize and output ice table-related variables and in the launching script.
JBC

File size: 6.0 KB
Line 
1MODULE pemredem
2
3!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
4!!!
5!!! Purpose: Write specific netcdf restart for the PEM
6!!!
7!!!
8!!! Author: LL, inspired by phyredem from the PCM
9!!!
10!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
11
12implicit none
13
14!=======================================================================
15contains
16!=======================================================================
17
18SUBROUTINE pemdem0(filename,lonfi,latfi,cell_area,ngrid,nslope,def_slope,subslope_dist)
19
20! create physics restart file and write time-independent variables
21use comsoil_h_PEM, only: mlayer_PEM
22use iostart_PEM,   only: open_restartphy, close_restartphy, put_var, put_field, length
23
24implicit none
25
26character(*),                  intent(in) :: filename
27integer,                       intent(in) :: ngrid, nslope
28real, dimension(ngrid),        intent(in) :: lonfi, latfi
29real, dimension(ngrid),        intent(in) :: cell_area     ! boundaries for bining of the slopes
30real, dimension(nslope + 1),   intent(in) :: def_slope     ! boundaries for bining of the slopes
31real, dimension(ngrid,nslope), intent(in) :: subslope_dist ! undermesh statistics
32
33! Create physics start file
34call open_restartphy(filename)
35
36! Write the mid-layer depths
37call put_var("soildepth","Soil mid-layer depth",mlayer_PEM)
38
39! Write longitudes
40call put_field("longitude","Longitudes of physics grid",lonfi)
41
42! Write latitudes
43call put_field("latitude","Latitudes of physics grid",latfi)
44
45! Write mesh areas
46call put_field("area","Mesh area",cell_area)
47
48! Multidimensionnal variables (nopcm undermesh slope statistics)
49call put_var("def_slope","slope criterium stages",def_slope)
50call put_field("subslope_dist","under mesh slope distribution",subslope_dist)
51
52! Close file
53call close_restartphy
54
55END SUBROUTINE pemdem0
56
57!=======================================================================
58
59SUBROUTINE pemdem1(filename,i_myear,nsoil_PEM,ngrid,nslope,tsoil_slope_PEM,inertiesoil_slope_PEM, &
60                   icetable_depth,icetable_thickness,ice_porefilling,m_co2_regolith,m_h2o_regolith,h2o_ice,stratif)
61
62! write time-dependent variable to restart file
63use iostart_PEM,   only: open_restartphy, close_restartphy, put_var, put_field
64use comsoil_h_PEM, only: inertiedat_PEM, soil_pem
65use time_evol_mod, only: year_bp_ini, convert_years
66use layering_mod,  only: layering, nb_str_max, stratif2array, print_layering, layering_algo
67
68implicit none
69
70#ifndef CPP_STD
71    include "callkeys.h"
72#endif
73
74character(*),                            intent(in) :: filename
75integer,                                 intent(in) :: nsoil_PEM, ngrid, nslope
76real,                                    intent(in) :: i_myear
77real, dimension(ngrid,nsoil_PEM,nslope), intent(in) :: tsoil_slope_PEM       ! under mesh bining according to slope
78real, dimension(ngrid,nsoil_PEM,nslope), intent(in) :: inertiesoil_slope_PEM ! under mesh bining according to slope
79real, dimension(ngrid,nslope),           intent(in) :: icetable_depth       ! under mesh bining according to slope
80real, dimension(ngrid,nslope),           intent(in) :: icetable_thickness   ! under mesh bining according to slope
81real, dimension(ngrid,nsoil_PEM,nslope), intent(in) :: ice_porefilling       ! under mesh bining according to slope
82real, dimension(ngrid,nsoil_PEM,nslope), intent(in) :: m_co2_regolith, m_h2o_regolith
83real, dimension(ngrid,nslope),           intent(in) :: h2o_ice
84type(layering), dimension(ngrid,nslope), intent(in) :: stratif ! Stratification (layerings)
85
86integer                               :: islope
87character(2)                          :: num
88real                                  :: Year          ! Year of the simulation
89real, dimension(:,:,:,:), allocatable :: stratif_array ! Array for stratification (layerings)
90
91! Open file
92call open_restartphy(filename)
93
94! First variable to write must be "Time", in order to correctly
95! set time counter in file
96Year = (year_bp_ini + i_myear)*convert_years
97call put_var("Time","Year of simulation",Year)
98call put_field('h2o_ice','h2o_ice',h2o_ice,Year)
99
100if (layering_algo) then
101    allocate(stratif_array(ngrid,nslope,nb_str_max,6))
102    call stratif2array(stratif,ngrid,nslope,stratif_array)
103    do islope = 1,nslope
104        write(num,fmt='(i2.2)') islope
105        call put_field('stratif_slope'//num//'_thickness','Layering thickness',stratif_array(:,islope,:,1),Year)
106        call put_field('stratif_slope'//num//'_top_elevation','Layering top elevation',stratif_array(:,islope,:,2),Year)
107        call put_field('stratif_slope'//num//'_co2ice_volfrac','Layering CO2 ice volume fraction',stratif_array(:,islope,:,3),Year)
108        call put_field('stratif_slope'//num//'_h2oice_volfrac','Layering H2O ice volume fraction',stratif_array(:,islope,:,4),Year)
109        call put_field('stratif_slope'//num//'_dust_volfrac','Layering dust volume fraction',stratif_array(:,islope,:,5),Year)
110        call put_field('stratif_slope'//num//'_air_volfrac','Layering air volume fraction',stratif_array(:,islope,:,6),Year)
111    enddo
112    deallocate(stratif_array)
113endif
114
115if (soil_pem) then
116  ! Multidimensionnal variables (undermesh slope statistics)
117    do islope = 1,nslope
118        write(num,fmt='(i2.2)') islope
119        call put_field("tsoil_PEM_slope"//num,"Soil temperature by slope type",tsoil_slope_PEM(:,:,islope),Year)
120        call put_field("TI_PEM_slope"//num,"Soil Thermal Inertia by slope type",inertiesoil_slope_PEM(:,:,islope),Year)
121        call put_field("mco2_reg_ads_slope"//num, "Mass of co2 adsorbded in the regolith",m_co2_regolith(:,:,islope),Year)
122        call put_field("mh2o_reg_ads_slope"//num, "Mass of h2o adsorbded in the regolith",m_h2o_regolith(:,:,islope),Year)
123    enddo
124    call put_field("icetable_depth","Depth of ice table",icetable_depth,Year)
125    call put_field("icetable_thickness","Depth of ice table",icetable_thickness,Year)
126    call put_field("ice_porefilling","Subsurface ice pore filling",ice_porefilling,Year)
127    call put_field("inertiedat_PEM","Thermal inertie of PEM ",inertiedat_PEM,Year)
128endif ! soil_pem
129
130! Close file
131call close_restartphy
132
133END SUBROUTINE pemdem1
134
135END MODULE pemredem
Note: See TracBrowser for help on using the repository browser.