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

Last change on this file since 3203 was 3172, checked in by llange, 18 months ago

PEM
Fixing bug when rewritting the startfi.nc for the PCM: fluxgeo, read in the run_PEM.def, was not written in the startfi.

LL

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