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

Last change on this file since 3241 was 3206, checked in by jbclement, 17 months ago

PEM:
Cleanings of unused variables/arguments and bad type conversions.
JBC

File size: 4.4 KB
RevLine 
[3149]1MODULE pemredem
2
[2855]3!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
4!!!
[3149]5!!! Purpose: Write specific netcdf restart for the PEM
[2855]6!!!
[3149]7!!!
8!!! Author: LL, inspired by phyredem from the PCM
9!!!
[2855]10!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
[2794]11
12implicit none
13
[3149]14!=======================================================================
[2794]15contains
[3149]16!=======================================================================
[2794]17
[3206]18SUBROUTINE pemdem0(filename,lonfi,latfi,cell_area,ngrid,nslope,def_slope,subslope_dist)
[3149]19
[2794]20! create physics restart file and write time-independent variables
[3206]21use comsoil_h_PEM, only: mlayer_PEM
22use iostart_PEM,   only: open_restartphy, close_restartphy, put_var, put_field, length
[3149]23
24implicit none
[2794]25
[3149]26character(*),                  intent(in) :: filename
[3206]27integer,                       intent(in) :: ngrid, nslope
[3149]28real, dimension(ngrid),        intent(in) :: lonfi, latfi
29real, dimension(ngrid),        intent(in) :: cell_area     ! boundaries for bining of the slopes
30real, dimension(ngrid + 1),    intent(in) :: def_slope     ! boundaries for bining of the slopes
31real, dimension(ngrid,nslope), intent(in) :: subslope_dist ! undermesh statistics
[2794]32
[3149]33! Create physics start file
34call open_restartphy(filename)
[2794]35
[3149]36! Write the mid-layer depths
37call put_var("soildepth","Soil mid-layer depth",mlayer_PEM)
[2888]38
[3149]39! Write longitudes
40call put_field("longitude","Longitudes of physics grid",lonfi)
[2888]41
[3149]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                   ice_table_depth,ice_table_thickness,m_co2_regolith,m_h2o_regolith,h2o_ice)
61
62! write time-dependent variable to restart file
63use iostart_PEM,   only: open_restartphy, close_restartphy, put_var, put_field
[3206]64use comsoil_h_PEM, only: inertiedat_PEM, soil_pem
[3149]65use time_evol_mod, only: year_bp_ini, convert_years
66
67implicit none
68
[2842]69#ifndef CPP_STD
[3149]70    include "callkeys.h"
[2842]71#endif
[2794]72
[3149]73character(*),                            intent(in) :: filename
74integer,                                 intent(in) :: nsoil_PEM, ngrid, nslope, i_myear
75real, dimension(ngrid,nsoil_PEM,nslope), intent(in) :: tsoil_slope_PEM       ! under mesh bining according to slope
76real, dimension(ngrid,nsoil_PEM,nslope), intent(in) :: inertiesoil_slope_PEM ! under mesh bining according to slope
77real, dimension(ngrid,nslope),           intent(in) :: ice_table_depth       ! under mesh bining according to slope
78real, dimension(ngrid,nslope),           intent(in) :: ice_table_thickness   ! under mesh bining according to slope
79real, dimension(ngrid,nsoil_PEM,nslope), intent(in) :: m_co2_regolith, m_h2o_regolith
80real, dimension(ngrid,nslope),           intent(in) :: h2o_ice
[2794]81
[3149]82integer       :: islope
83character(2)  :: num
84real          :: Year ! Year of the simulation
[3039]85
[3149]86! Open file
87call open_restartphy(filename)
[2794]88
[3149]89! First variable to write must be "Time", in order to correctly
90! set time counter in file
91Year = (year_bp_ini + i_myear)*convert_years
92call put_var("Time","Year of simulation",Year)
93call put_field('h2o_ice','h2o_ice',h2o_ice,Year)
[2794]94
[3149]95if (soil_pem) then
96  ! Multidimensionnal variables (undermesh slope statistics)
97    do islope = 1,nslope
98        write(num,fmt='(i2.2)') islope
99        call put_field("tsoil_PEM_slope"//num,"Soil temperature by slope type",tsoil_slope_PEM(:,:,islope),Year)
100        call put_field("TI_PEM_slope"//num,"Soil Thermal Inertia by slope type",inertiesoil_slope_PEM(:,:,islope),Year)
101        call put_field("mco2_reg_ads_slope"//num, "Mass of co2 adsorbded in the regolith",m_co2_regolith(:,:,islope),Year)
102        call put_field("mh2o_reg_ads_slope"//num, "Mass of h2o adsorbded in the regolith",m_h2o_regolith(:,:,islope),Year)
103    enddo
104    call put_field("ice_table_depth","Depth of ice table",ice_table_depth,Year)
105    call put_field("ice_table_thickness","Depth of ice table",ice_table_thickness,Year)
106    call put_field("inertiedat_PEM","Thermal inertie of PEM ",inertiedat_PEM,Year)
107endif ! soil_pem
[2794]108
[3149]109! Close file
110call close_restartphy
[2794]111
[3149]112END SUBROUTINE pemdem1
[2961]113
[3149]114END MODULE pemredem
Note: See TracBrowser for help on using the repository browser.