source: trunk/LMDZ.COMMON/libf/evolution/pem.F90 @ 3939

Last change on this file since 3939 was 3938, checked in by jbclement, 3 months ago

PEM:

  • Correction of the process to balance the H2O flux from and into the atmosphere accross reservoirs: (i) computation of H2O amount going from/in the atmosphere is corrected (missing dt and wrong condition); (ii) computation of the scaling factor is corrected because we need to balance all the icetable/adsorbed/surface ice flux but it affects only sublimating/condensing surface ice flux; (iii) process of balancing is modified to be applied to every flux by scaling instead of applying it only to positive or negative flux by trimming; (iv) balance is now fully processed by the tendency before evolving the ice. This is done through dedicated subroutines.
  • Correction of ice conversion between the layering algorithm and PEM ice variables (density factor was missing).
  • Addition of H2O flux balance and related stopping criterion for the layering algorithm.
  • Few updates for files in the deftank to be in line with the wiki Planeto pages.

JBC

File size: 67.4 KB
RevLine 
[2779]1!------------------------
[3028]2! I   Initialization
[3161]3!    I_a Read the "run.def"
[3317]4!    I_b Read the "start.nc" and "startfi.nc"
[2835]5!    I_c Subslope parametrisation
[3161]6!    I_d Read the PCM data and convert them to the physical grid
[3028]7!    I_e Initialization of the PEM variable and soil
[3143]8!    I_f Compute tendencies
[3571]9!    I_g Compute global surface pressure
[3384]10!    I_h Read the "startpem.nc"
[2835]11!    I_i Compute orbit criterion
[2779]12
13! II  Run
[3028]14!    II_a Update pressure, ice and tracers
[3149]15!    II_b Evolution of ice
16!    II_c Flow of glaciers
[2835]17!    II_d Update surface and soil temperatures
[3088]18!    II_e Outputs
19!    II_f Update the tendencies
20!    II_g Checking the stopping criterion
[2779]21
22! III Output
[3770]23!    III_a Update surface values for the PCM start files
[3317]24!    III_b Write the "restart.nc" and "restartfi.nc"
[3161]25!    III_c Write the "restartpem.nc"
[2779]26!------------------------
27
28PROGRAM pem
29
[3161]30use phyetat0_mod,               only: phyetat0
31use phyredem,                   only: physdem0, physdem1
32use netcdf,                     only: nf90_open, NF90_NOWRITE, nf90_get_var, nf90_inq_varid, nf90_close
33use turb_mod,                   only: q2, wstar
[3206]34use comslope_mod,               only: nslope, def_slope, def_slope_mean, subslope_dist, iflat, ini_comslope_h
[3161]35use logic_mod,                  only: iflag_phys
36use mod_const_mpi,              only: COMM_LMDZ
[3028]37use infotrac
[3161]38use geometry_mod,               only: latitude_deg
39use conf_pem_mod,               only: conf_pem
40use pemredem,                   only: pemdem0, pemdem1
41use glaciers_mod,               only: flow_co2glaciers, flow_h2oglaciers, co2ice_flow, h2oice_flow, inf_h2oice_threshold, &
[3571]42                                      metam_h2oice_threshold, metam_co2ice_threshold, metam_h2oice, metam_co2ice, computeTcondCO2
[3938]43use stopping_crit_mod,          only: stopping_crit_h2o_ice, stopping_crit_co2, stopping_crit_h2o
[3161]44use constants_marspem_mod,      only: alpha_clap_co2, beta_clap_co2, alpha_clap_h2o, beta_clap_h2o, m_co2, m_noco2
[3938]45use evol_ice_mod,               only: evol_co2_ice, evol_h2o_ice, balance_h2oice_reservoirs
[3571]46use comsoil_h_PEM,              only: soil_pem, ini_comsoil_h_PEM, end_comsoil_h_PEM, nsoilmx_PEM, &
47                                      TI_PEM,               & ! Soil thermal inertia
[3330]48                                      tsoil_PEM, layer_PEM, & ! Soil temp, number of subsurface layers, soil mid layer depths
49                                      fluxgeo                 ! Geothermal flux for the PEM and PCM
[3161]50use adsorption_mod,             only: regolith_adsorption, adsorption_pem,        & ! Bool to check if adsorption, main subroutine
51                                      ini_adsorption_h_PEM, end_adsorption_h_PEM, & ! Allocate arrays
[3571]52                                      co2_adsorbed_phys, h2o_adsorbed_phys          ! Mass of co2 and h2O adsorbed
[3498]53use time_evol_mod,              only: dt, evol_orbit_pem, Max_iter_pem, convert_years, year_bp_ini
[3161]54use orbit_param_criterion_mod,  only: orbit_param_criterion
55use recomp_orb_param_mod,       only: recomp_orb_param
[3493]56use ice_table_mod,              only: icetable_depth, icetable_thickness, end_ice_table, ice_porefilling, &
57                                      ini_ice_table, icetable_equilibrium, icetable_dynamic, computeice_table_equilibrium, compute_massh2o_exchange_ssi
[3161]58use soil_thermalproperties_mod, only: update_soil_thermalproperties
[3206]59use time_phylmdz_mod,           only: daysec, dtphys
[3161]60use abort_pem_mod,              only: abort_pem
61use soil_settings_PEM_mod,      only: soil_settings_PEM
62use compute_tend_mod,           only: compute_tend
63use info_PEM_mod,               only: info_PEM
[3571]64use get_timelen_PCM_mod,        only: get_timelen_PCM
[3161]65use pemetat0_mod,               only: pemetat0
66use read_data_PCM_mod,          only: read_data_PCM
[3610]67use recomp_tend_mod,            only: recomp_tend_co2, recomp_tend_h2o
[3553]68use compute_soiltemp_mod,       only: compute_tsoil_pem, shift_tsoil2surf
[3181]69use writediagpem_mod,           only: writediagpem, writediagsoilpem
[3207]70use co2condens_mod,             only: CO2cond_ps
[3850]71use layering_mod,               only: layering, del_layering, make_layering, layering_algo, subsurface_ice_layering, &
[3938]72                                      rho_co2ice, rho_h2oice, ptrarray,                                              &
73                                      stratum, get_nb_str_max, nb_str_max, is_dust_lag, is_co2ice_str, is_h2oice_str
[3527]74use dyn_ss_ice_m_mod,           only: dyn_ss_ice_m
[3836]75use parse_args_mod,             only: parse_args
76use job_timelimit_mod,          only: timelimit, antetime, timewall
[3638]77use paleoclimate_mod,           only: h2o_ice_depth, zdqsdif_ssi_tot
[3907]78use surf_temp,                  only: update_tsurf_nearest_baresoil
[2985]79
[2842]80#ifndef CPP_STD
[3206]81    use comsoil_h,          only: tsoil, nsoilmx, ini_comsoil_h, inertiedat, mlayer, inertiesoil, flux_geo, nqsoil, qsoil
[3584]82    use surfdat_h,          only: tsurf, qsurf, emis, emissiv, emisice, ini_surfdat_h,   &
83                                  albedodat, albedice, albedo_h2o_frost, albedo_h2o_cap, &
84                                  zmea, zstd, zsig, zgam, zthe, frost_albedo_threshold,  &
85                                  watercap, watercaptag, perennial_co2ice, albedo_perennialco2
[3028]86    use dimradmars_mod,     only: totcloudfrac, albedo
87    use dust_param_mod,     only: tauscaling
[3149]88    use tracer_mod,         only: noms, igcm_h2o_ice, igcm_co2, mmol, igcm_h2o_vap ! Tracer names and molar masses
[3028]89    use mod_phys_lmdz_para, only: is_parallel, is_sequential, is_mpi_root, is_omp_root, is_master
[3096]90    use planete_h,          only: aphelie, periheli, year_day, peri_day, obliquit, iniorbit
[3143]91    use surfini_mod,        only: surfini
[3618]92    use comcstfi_h,         only: mugaz
[2842]93#else
[3028]94    use tracer_h,           only: noms, igcm_h2o_ice, igcm_co2 ! Tracer names
95    use phys_state_var_mod, only: cloudfrac, totcloudfrac, albedo_snow_SPECTV,HICE,RNAT,   &
96                                  PCTSRF_SIC, TSLAB, TSEA_ICE, SEA_ICE, ALBEDO_BAREGROUND, &
97                                  ALBEDO_CO2_ICE_SPECTV, phys_state_var_init
[3039]98    use aerosol_mod,        only: iniaerosol
99    use planete_mod,        only: apoastr, periastr, year_day, peri_day, obliquit
[3779]100    use comcstfi_mod,       only: pi, rad, g, r, cpp, rcp, mugaz
[2842]101#endif
[2985]102
[3028]103#ifndef CPP_1D
[3779]104    use comconst_mod,             only: pi, rad, g, r, cpp, rcp => kappa
[3076]105    use iniphysiq_mod,            only: iniphysiq
106    use control_mod,              only: iphysiq, day_step, nsplit_phys
[3019]107#else
[3779]108    use comcstfi_h,               only: pi, rad, g, r, cpp, rcp
[3386]109    use time_phylmdz_mod,         only: iphysiq, steps_per_sol
[3028]110    use regular_lonlat_mod,       only: init_regular_lonlat
111    use physics_distribution_mod, only: init_physics_distribution
112    use mod_grid_phy_lmdz,        only: regular_lonlat
[3065]113    use init_testphys1d_mod,      only: init_testphys1d
114    use comvert_mod,              only: ap, bp
[3076]115    use writerestart1D_mod,       only: writerestart1D
[2980]116#endif
[2835]117
[3076]118implicit none
[2980]119
[3028]120include "dimensions.h"
121include "paramet.h"
122include "comgeom.h"
123include "iniprint.h"
[2779]124
[3096]125! Same variable names as in the PCM
[3770]126integer, parameter :: ngridmx = 2 + (jjm - 1)*iim - 1/jjm
[3065]127integer, parameter :: nlayer = llm ! Number of vertical layer
128integer            :: ngrid        ! Number of physical grid points
129integer            :: nq           ! Number of tracer
130integer            :: day_ini      ! First day of the simulation
131real               :: pday         ! Physical day
[3149]132real               :: time_phys    ! Same as in PCM
133real               :: ptimestep    ! Same as in PCM
[3571]134real               :: ztime_fin    ! Same as in PCM
[2794]135
[3571]136! Variables to read "start.nc"
[3317]137character(*), parameter :: start_name = "start.nc" ! Name of the file used to initialize the PEM
[2779]138
[3028]139! Dynamic variables
[3065]140real, dimension(ip1jm,llm)          :: vcov          ! vents covariants
141real, dimension(ip1jmp1,llm)        :: ucov          ! vents covariants
[3571]142real, dimension(ip1jmp1,llm)        :: teta          ! Potential temperature
[3065]143real, dimension(:,:,:), allocatable :: q             ! champs advectes
[3603]144real, dimension(:),     allocatable :: pdyn          ! pressure for the dynamic grid
[3571]145real, dimension(:),     allocatable :: ps_start      ! surface pressure in the start file
146real, dimension(:),     allocatable :: ps_start0     ! surface pressure in the start file at the beginning
[3673]147real, dimension(:),     allocatable :: ps_avg        ! (ngrid) Average surface pressure
[3571]148real, dimension(:),     allocatable :: ps_dev        ! (ngrid x timelen) Surface pressure deviation
149real, dimension(:,:),   allocatable :: ps_timeseries ! (ngrid x timelen) Instantaneous surface pressure
150real, dimension(ip1jmp1,llm)        :: masse         ! Air mass
[3065]151real, dimension(ip1jmp1)            :: phis          ! geopotentiel au sol
[3028]152real                                :: time_0
[2779]153
[3028]154! Variables to read starfi.nc
[3317]155character(*), parameter :: startfi_name = "startfi.nc" ! Name of the file used to initialize the PEM
[3143]156character(2)            :: str2
[3206]157integer                 :: ncid, status                           ! Variable for handling opening of files
158integer                 :: lonvarid, latvarid, areavarid, sdvarid ! Variable ID for Netcdf files
159integer                 :: apvarid, bpvarid                       ! Variable ID for Netcdf files
[2794]160
[3028]161! Variables to read starfi.nc and write restartfi.nc
[3143]162real, dimension(:), allocatable :: longitude     ! Longitude read in startfi_name and written in restartfi
163real, dimension(:), allocatable :: latitude      ! Latitude read in startfi_name and written in restartfi
164real, dimension(:), allocatable :: cell_area     ! Cell_area read in startfi_name and written in restartfi
[3571]165real                            :: total_surface ! Total surface of the planet
[2897]166
[3778]167! Variables for h2o ice evolution
[3327]168real, dimension(:,:),    allocatable  :: h2o_ice              ! h2o ice in the PEM
[3571]169real, dimension(:,:),    allocatable  :: d_h2oice             ! physical point x slope field: Tendency of evolution of perennial h2o ice
[3327]170real, dimension(:,:,:),  allocatable  :: min_h2o_ice          ! Minima of h2o ice at each point for the PCM years [kg/m^2]
171real                                  :: h2oice_ini_surf      ! Initial surface of sublimating h2o ice
[3571]172logical, dimension(:,:), allocatable  :: is_h2oice_sublim_ini ! Logical array to know if h2o ice is sublimating
173real                                  :: ps_avg_global_ini    ! constant: Global average pressure at initialization [Pa]
174real                                  :: ps_avg_global_old    ! constant: Global average pressure of previous time step
175real                                  :: ps_avg_global_new    ! constant: Global average pressure of current time step
176real,   dimension(:,:),   allocatable :: zplev_new            ! Grid points x Atmospheric field: mass of the atmospheric layers in the pem at current time step [kg/m^2]
177real,   dimension(:,:),   allocatable :: zplev_start0         ! Grid points x Atmospheric field: mass of the atmospheric layers in the start [kg/m^2]
178real,   dimension(:,:,:), allocatable :: zplev_timeseries_new ! Grid points x Atmospheric x Time: same as zplev_new, but in times series [kg/m ^2]
179real,   dimension(:,:,:), allocatable :: zplev_timeseries_old ! same but with the time series, for previous time step
[3327]180integer                               :: stopPEM              ! which criterion is reached? 0 = no stopping; 1 = h2o ice surf; 2 = no h2o ice; 3 = co2 ice surf; 4 = ps; 5 = orb param; 6 = end of simu
[3571]181real                                  :: A, B, mmean          ! Molar mass: intermediate A, B for computations of the  mean molar mass of the layer [mol/kg]
182real,   dimension(:,:),   allocatable :: q_h2o_PEM_phys       ! Grid points x Times: h2o mass mixing ratio computed in the PEM, first value comes from PCM [kg/kg]
[3327]183integer                               :: timelen              ! # time samples
184real                                  :: extra_mass           ! Intermediate variables Extra mass of a tracer if it is greater than 1
[3938]185real, dimension(:,:),    allocatable  :: d_h2oice_new         ! Adjusted tendencies to keep balance between donor and recipient
186real                                  :: S_atm_2_h2o, S_h2o_2_atm, S_atm_2_h2oice, S_h2oice_2_atm ! Variables to balance H2O
[2779]187
[3778]188! Variables for co2 ice evolution
[3571]189real,    dimension(:,:), allocatable :: co2_ice                ! co2 ice in the PEM
190real,    dimension(:,:), allocatable :: d_co2ice               ! physical point x slope field: Tendency of evolution of perennial co2 ice over a year
191real,    dimension(:,:), allocatable :: d_co2ice_ini           ! physical point x slope field: Tendency of evolution of perennial co2 ice over a year in the PCM
192logical, dimension(:,:), allocatable :: is_co2ice_ini          ! Was there co2 ice initially in the PEM?
193real,  dimension(:,:,:), allocatable :: min_co2_ice            ! Minimum of co2 ice at each point for the first year [kg/m^2]
194real                                 :: co2ice_sublim_surf_ini ! Initial surface of sublimating co2 ice
195logical, dimension(:,:), allocatable :: is_co2ice_sublim_ini   ! Logical array to know if co2 ice is sublimating
196real,    dimension(:,:), allocatable :: vmr_co2_PCM            ! Grid points x Times co2 volume mixing ratio retrieve from the PCM [m^3/m^3]
197real,    dimension(:,:), allocatable :: vmr_co2_PEM_phys       ! Grid points x Times co2 volume mixing ratio used in the PEM
198real,    dimension(:,:), allocatable :: q_co2_PEM_phys         ! Grid points x Times co2 mass mixing ratio in the first layer computed in the PEM, first value comes from PCM [kg/kg]
[3779]199real(kind = 16)                      :: totmass_co2ice, totmass_atmco2         ! Current total CO2 masses
200real(kind = 16)                      :: totmass_co2ice_ini, totmass_atmco2_ini ! Initial total CO2 masses
[3130]201
[3778]202! Variables for the evolution of layered layerings_map
203type(layering), dimension(:,:), allocatable :: layerings_map     ! Layering for each grid point and slope
[3770]204type(ptrarray), dimension(:,:), allocatable :: current           ! Current active stratum in the layering
205logical,        dimension(:,:), allocatable :: new_str, new_lag  ! Flags for the layering algorithm
[3784]206real,           dimension(:,:), allocatable :: h2o_ice_depth_old ! Old depth of subsurface ice layer
[3297]207
[3028]208! Variables for slopes
[3591]209integer(kind = 1), dimension(:,:), allocatable :: flag_co2flow ! (ngrid,nslope): Flag where there is a CO2 glacier flow
210integer(kind = 1), dimension(:,:), allocatable :: flag_h2oflow ! (ngrid,nslope): Flag where there is a H2O glacier flow
[2779]211
[3028]212! Variables for surface and soil
[3779]213real, dimension(:,:),     allocatable :: tsurf_avg                          ! Grid points x Slope field: Average surface temperature [K]
214real, dimension(:,:),     allocatable :: tsurf_dev                          ! Grid points x Slope field: Surface temperature deviation [K]
215real, dimension(:,:),     allocatable :: tsurf_avg_yr1                      ! Grid points x Slope field: Average surface temperature of first call of the PCM [K]
216real, dimension(:,:,:),   allocatable :: tsoil_avg                          ! Grid points x Soil x Slope field: Average Soil Temperature [K]
217real, dimension(:,:),     allocatable :: tsoil_avg_old                      ! Grid points x Soil field: Average Soil Temperature at the previous time step [K]
218real, dimension(:,:,:),   allocatable :: tsoil_dev                          ! Grid points x Soil x Slope field: Soil temperature deviation [K]
219real, dimension(:,:,:,:), allocatable :: tsoil_timeseries                   ! Grid points x Soil x Slope x Times field: Soil temperature timeseries [K]
220real, dimension(:,:,:,:), allocatable :: tsoil_PEM_timeseries               ! Grid points x Soil x Slope x Times field: Soil temperature timeseries for PEM [K]
221real, dimension(:,:,:,:), allocatable :: tsoil_PEM_timeseries_old           ! Grid points x Soil x Slope x Times field: Soil temperature timeseries for PEM at the previous time step [K]
222real, dimension(:,:,:,:), allocatable :: watersoil_density_timeseries       ! Grid points x Soil x Slope x Times Water soil density timeseries [kg /m^3]
223real, dimension(:,:),     allocatable :: watersurf_density_avg              ! Grid points x Slope: Average water surface density [kg/m^3]
224real, dimension(:,:,:,:), allocatable :: watersoil_density_PEM_timeseries   ! Grid points x Soil x Slope x Times: Water soil density timeseries for PEM [kg/m^3]
225real, dimension(:,:,:),   allocatable :: watersoil_density_PEM_avg          ! Grid points x Soil x Slopes: Average water soil density [kg/m^3]
226real, dimension(:),       allocatable :: delta_co2_adsorbed                 ! Physics: quantity of CO2 that is exchanged because of adsorption / desorption [kg/m^2]
227real, dimension(:),       allocatable :: delta_h2o_adsorbed                 ! Physics: quantity of H2O that is exchanged because of adsorption / desorption [kg/m^2]
228real(kind = 16)                       :: totmass_adsco2, totmass_adsco2_ini ! Total mass of CO2 that is exchanged because of adsorption / desoprtion over the planets [kg]
229real                                  :: totmass_adsh2o                     ! Total mass of H2O that is exchanged because of adsorption / desoprtion over the planets [kg]
230logical, dimension(:,:),  allocatable :: co2ice_disappeared                 ! logical to check if a co2 ice reservoir already disappeared at a previous timestep
231real, dimension(:,:),     allocatable :: icetable_thickness_old             ! ngrid x nslope: Thickness of the ice table at the previous iteration [m]
232real, dimension(:,:,:),   allocatable :: ice_porefilling_old                ! ngrid x nslope: Ice pore filling at the previous iteration [m]
233real, dimension(:),       allocatable :: delta_h2o_icetablesublim           ! ngrid x Total mass of the H2O that has sublimated / condenses from the ice table [kg]
234real, dimension(:),       allocatable :: porefill                           ! Pore filling (output) to compute the dynamic ice table
235real                                  :: ssi_depth                          ! Ice table depth (output) to compute the dynamic ice table
236real, dimension(:,:),     allocatable :: zshift_surf                        ! Elevation shift for the surface [m]
237real, dimension(:,:),     allocatable :: zlag                               ! Newly built lag thickness [m]
238real, dimension(:,:),     allocatable :: icetable_depth_old                 ! Old depth of the ice table
[3065]239
[3028]240! Some variables for the PEM run
[3616]241real, parameter       :: year_step = 1   ! Timestep for the pem
242real                  :: i_myear_leg     ! Number of iteration
243real                  :: n_myear_leg     ! Maximum number of iterations before stopping
244real                  :: i_myear         ! Global number of Martian years of the chained simulations
245real                  :: n_myear         ! Maximum number of Martian years of the chained simulations
246real                  :: timestep        ! Timestep [s]
247integer(kind = 8)     :: cr              ! Number of clock ticks per second (count rate)
248integer(kind = 8)     :: c1, c2          ! Counts of processor clock
249character(8)          :: date
250character(10)         :: time
251character(5)          :: zone
252integer, dimension(8) :: values
253character(128)        :: dir = ' '
254character(32)         :: logname = ' '
255character(32)         :: hostname = ' '
[2779]256
[2842]257#ifdef CPP_STD
[3571]258    real                                :: frost_albedo_threshold = 0.05 ! Frost albedo threeshold to convert fresh frost to old ice
259    real                                :: albedo_h2o_frost              ! Albedo of h2o frost
[3143]260    real, dimension(:),     allocatable :: tsurf_read_generic            ! Temporary variable to do the subslope transfert dimension when reading form generic
261    real, dimension(:,:),   allocatable :: qsurf_read_generic            ! Temporary variable to do the subslope transfert dimension when reading form generic
262    real, dimension(:,:),   allocatable :: tsoil_read_generic            ! Temporary variable to do the subslope transfert dimension when reading form generic
263    real, dimension(:),     allocatable :: emis_read_generic             ! Temporary variable to do the subslope transfert dimension when reading form generic
264    real, dimension(:,:),   allocatable :: albedo_read_generic           ! Temporary variable to do the subslope transfert dimension when reading form generic
[3065]265    real, dimension(:,:),   allocatable :: tsurf                         ! Subslope variable, only needed in the GENERIC case
266    real, dimension(:,:,:), allocatable :: qsurf                         ! Subslope variable, only needed in the GENERIC case
267    real, dimension(:,:,:), allocatable :: tsoil                         ! Subslope variable, only needed in the GENERIC case
268    real, dimension(:,:),   allocatable :: emis                          ! Subslope variable, only needed in the GENERIC case
269    real, dimension(:,:),   allocatable :: watercap                      ! Subslope variable, only needed in the GENERIC case =0 no watercap in generic model
[3068]270    logical, dimension(:),  allocatable :: watercaptag                   ! Subslope variable, only needed in the GENERIC case =false no watercaptag in generic model
[3065]271    real, dimension(:,:,:), allocatable :: albedo                        ! Subslope variable, only needed in the GENERIC case
272    real, dimension(:,:,:), allocatable :: inertiesoil                   ! Subslope variable, only needed in the GENERIC case
[2842]273#endif
274
[2980]275#ifdef CPP_1D
[3143]276    integer            :: nsplit_phys
277    integer, parameter :: jjm_value = jjm - 1
[3386]278    integer            :: day_step
[3065]279
280    ! Dummy variables to use the subroutine 'init_testphys1d'
[3912]281    logical                             :: therestart1D, therestartfi, ctrl_h2ovap, ctrl_h2oice
[3068]282    integer                             :: ndt, day0
[3912]283    real                                :: ptif, pks, day, gru, grv, relaxtime_h2ovap, relaxtime_h2oice
[3068]284    real, dimension(:),     allocatable :: zqsat
285    real, dimension(:,:,:), allocatable :: dq, dqdyn
[3912]286    real, dimension(nlayer)             :: play, w, qref_h2ovap, qref_h2oice
[3068]287    real, dimension(nlayer + 1)         :: plev
[2980]288#else
[3591]289    integer, parameter                :: jjm_value = jjm
290    real, dimension(:),   allocatable :: ap ! Coefficient ap read in start_name and written in restart
291    real, dimension(:),   allocatable :: bp ! Coefficient bp read in start_name and written in restart
[3628]292    real, dimension(:,:), allocatable :: p  ! Grid points x Atmosphere: pressure to recompute and write in restart (ip1jmp1,llmp1)
[2980]293#endif
294
[3028]295! Loop variables
[3584]296integer :: i, l, ig, nnq, t, islope, ig_loop, islope_loop, isoil, icap
[3851]297real    :: totmass_ini
[3599]298logical :: num_str
[2779]299
[3786]300! CODE
[3842]301! Elapsed time with system clock
302call system_clock(count_rate = cr)
303call system_clock(c1)
304
305! Parse command-line options
306call parse_args()
307
308! Header
[3649]309write(*,*) '  *    .          .   +     .    *        .  +      .    .       .      '
[3603]310write(*,*) '           +         _______  ________  ____    ____      *           + '
311write(*,*) ' +   .        *     |_   __ \|_   __  ||_   \  /   _|          .       *'
312write(*,*) '          .     .     | |__) | | |_ \_|  |   \/   |  *        *      .  '
313write(*,*) '       .              |  ___/  |  _| _   | |\  /| |      .        .     '
314write(*,*) '.  *          *      _| |_    _| |__/ | _| |_\/_| |_                  * '
[3649]315write(*,*) '            +       |_____|  |________||_____||_____|   +     .         '
316write(*,*) '  .      *          .   *       .   +       *          .        +      .'
[3603]317
[3616]318! Some user info
319call date_and_time(date,time,zone,values)
320call getcwd(dir)      ! Current directory
321call getlog(logname)  ! User name
322call hostnm(hostname) ! Machine/station name
323write(*,*)
[3673]324write(*,*) '********* PEM information *********'
[3616]325write(*,*) '+ User     : '//trim(logname)
326write(*,*) '+ Machine  : '//trim(hostname)
327write(*,*) '+ Directory: '//trim(dir)
328write(*,'(a,i2,a,i2,a,i4)') ' + Date     : ',values(3),'/',values(2),'/',values(1)
329write(*,'(a,i2,a,i2,a,i2,a)') ' + Time     : ',values(5),':',values(6),':',values(7)
330
[3028]331! Parallel variables
[2842]332#ifndef CPP_STD
[3028]333    is_sequential = .true.
334    is_parallel = .false.
335    is_mpi_root = .true.
336    is_omp_root = .true.
337    is_master = .true.
[2842]338#endif
[2779]339
[3065]340! Some constants
[3591]341day_ini = 0
342time_phys = 0.
[3028]343ngrid = ngridmx
[3591]344A = (1./m_co2 - 1./m_noco2)
345B = 1./m_noco2
[3028]346year_day = 669
347daysec = 88775.
[3426]348timestep = year_day*daysec*year_step
[2794]349
[3028]350!----------------------------- INITIALIZATION --------------------------
[2779]351!------------------------
[3028]352! I   Initialization
[3161]353!    I_a Read the "run.def"
[2779]354!------------------------
[3616]355write(*,*)
[3673]356write(*,*) '********* PEM initialization *********'
[3603]357write(*,*) '> Reading "run.def" (PEM)'
[2980]358#ifndef CPP_1D
[3383]359    dtphys = daysec/48. ! Dummy value (overwritten in phyetat0)
[3028]360    call conf_gcm(99,.true.)
361    call infotrac_init
362    nq = nqtot
363    allocate(q(ip1jmp1,llm,nqtot))
[3065]364    allocate(longitude(ngrid),latitude(ngrid),cell_area(ngrid))
[2980]365#else
[3603]366    allocate(q(1,llm,nqtot),pdyn(1))
[3065]367    allocate(longitude(1),latitude(1),cell_area(1))
[3129]368
[3143]369    therestart1D = .false. ! Default value
[3363]370    inquire(file = 'start1D.txt',exist = therestart1D)
[3129]371    if (.not. therestart1D) then
[3363]372        write(*,*) 'There is no "start1D.txt" file!'
[3129]373        error stop 'Initialization cannot be done for the 1D PEM.'
374    endif
[3143]375    therestartfi = .false. ! Default value
[3317]376    inquire(file = 'startfi.nc',exist = therestartfi)
[3129]377    if (.not. therestartfi) then
[3317]378        write(*,*) 'There is no "startfi.nc" file!'
[3129]379        error stop 'Initialization cannot be done for the 1D PEM.'
380    endif
381
[3603]382    write(*,*) '> Reading "start1D.txt" and "startfi.nc"'
[3609]383    call init_testphys1d('start1D.txt','startfi.nc',therestart1D,therestartfi,ngrid,nlayer,610.,nq,q,         &
384                         time_0,pdyn(1),ucov,vcov,teta,ndt,ptif,pks,dtphys,zqsat,dq,dqdyn,day0,day,gru,grv,w, &
[3912]385                         play,plev,latitude,longitude,cell_area,                                              &
386                         ctrl_h2ovap,relaxtime_h2ovap,qref_h2ovap,ctrl_h2oice,relaxtime_h2oice,qref_h2oice)
[3028]387    nsplit_phys = 1
[3399]388    day_step = steps_per_sol
[2980]389#endif
[2779]390
[3039]391call conf_pem(i_myear,n_myear)
[2779]392
[2835]393!------------------------
[3028]394! I   Initialization
[3571]395!    I_b Read of the "start.nc" and "starfi.nc"
[3028]396!------------------------
[3317]397! I_b.1 Read "start.nc"
[3603]398write(*,*) '> Reading "start.nc"'
[3591]399allocate(ps_start0(ngrid))
[2980]400#ifndef CPP_1D
[3603]401    allocate(pdyn(ip1jmp1))
402    call dynetat0(start_name,vcov,ucov,teta,q,masse,pdyn,phis,time_0)
[2779]403
[3603]404    call gr_dyn_fi(1,iip1,jjp1,ngridmx,pdyn,ps_start0)
[2897]405
[3589]406    call iniconst ! Initialization of dynamical constants (comconst_mod)
407    call inigeom ! Initialization of geometry
[2980]408
[3028]409    allocate(ap(nlayer + 1))
410    allocate(bp(nlayer + 1))
[3143]411    status = nf90_open(start_name,NF90_NOWRITE,ncid)
[3028]412    status = nf90_inq_varid(ncid,"ap",apvarid)
413    status = nf90_get_var(ncid,apvarid,ap)
414    status = nf90_inq_varid(ncid,"bp",bpvarid)
415    status = nf90_get_var(ncid,bpvarid,bp)
416    status = nf90_close(ncid)
[2779]417
[3591]418    ! Initialization of physics constants and variables (comcstfi_h)
[3571]419    call iniphysiq(iim,jjm,llm,(jjm - 1)*iim + 2,comm_lmdz,daysec,day_ini,dtphys/nsplit_phys,rlatu,rlatv,rlonu,rlonv,aire,cu,cv,rad,g,r,cpp,iflag_phys)
[2980]420#else
[3628]421    ps_start0(1) = pdyn(1)
[2980]422#endif
[3603]423deallocate(pdyn)
[2980]424
[3096]425! In the PCM, these values are given to the physic by the dynamic.
[3317]426! Here we simply read them in the "startfi.nc" file
[3571]427status = nf90_open(startfi_name,NF90_NOWRITE,ncid)
[3028]428status = nf90_inq_varid(ncid,"longitude",lonvarid)
429status = nf90_get_var(ncid,lonvarid,longitude)
430status = nf90_inq_varid(ncid,"latitude",latvarid)
431status = nf90_get_var(ncid,latvarid,latitude)
432status = nf90_inq_varid(ncid,"area",areavarid)
433status = nf90_get_var(ncid,areavarid,cell_area)
434status = nf90_inq_varid(ncid,"soildepth",sdvarid)
435status = nf90_get_var(ncid,sdvarid,mlayer)
436status = nf90_close(ncid)
[2963]437
[3317]438! I_b.2 Read the "startfi.nc"
[2779]439! First we read the initial state (starfi.nc)
[2842]440#ifndef CPP_STD
[3603]441    write(*,*) '> Reading "startfi.nc"'
[3143]442    call phyetat0(startfi_name,0,0,nsoilmx,ngrid,nlayer,nq,nqsoil,day_ini,time_phys,tsurf, &
[3149]443                  tsoil,albedo,emis,q2,qsurf,qsoil,tauscaling,totcloudfrac,wstar,          &
[3130]444                  watercap,perennial_co2ice,def_slope,def_slope_mean,subslope_dist)
[2779]445
[3070]446    ! Remove unphysical values of surface tracer
447    where (qsurf < 0.) qsurf = 0.
[2885]448
[3143]449    call surfini(ngrid,nslope,qsurf)
[2842]450#else
[3028]451    call phys_state_var_init(nq)
452    if (.not. allocated(noms)) allocate(noms(nq)) ! (because noms is an argument of physdem1 whether or not tracer is on)
453    call initracer(ngrid,nq)
454    call iniaerosol()
455    allocate(tsurf_read_generic(ngrid))
456    allocate(qsurf_read_generic(ngrid,nq))
457    allocate(tsoil_read_generic(ngrid,nsoilmx))
[3114]458    allocate(qsoil_read_generic(ngrid,nsoilmx,nqsoil,nslope))
[3028]459    allocate(emis_read_generic(ngrid))
[3591]460    allocate(albedo_read_generic(ngrid,2))
461    allocate(qsurf(ngrid,nq,1))
[3028]462    allocate(tsurf(ngrid,1))
463    allocate(tsoil(ngrid,nsoilmx,1))
464    allocate(emis(ngrid,1))
465    allocate(watercap(ngrid,1))
466    allocate(watercaptag(ngrid))
467    allocate(albedo(ngrid,2,1))
468    allocate(inertiesoil(ngrid,nsoilmx,1))
[3143]469    call phyetat0(.true.,ngrid,nlayer,startfi_name,0,0,nsoilmx,nq,nqsoil,day_ini,time_phys, &
[3149]470                  tsurf_read_generic,tsoil_read_generic,emis_read_generic,q2,               &
471                  qsurf_read_generic,qsoil_read_generic,cloudfrac,totcloudfrac,hice,        &
[3114]472                  rnat,pctsrf_sic,tslab,tsea_ice,sea_ice)
[3065]473    call surfini(ngrid,nq,qsurf_read_generic,albedo_read_generic,albedo_bareground,albedo_snow_SPECTV,albedo_co2_ice_SPECTV)
[2842]474
[3028]475    nslope = 1
476    call ini_comslope_h(ngrid,1)
[2842]477
[3149]478    qsurf(:,:,1) = qsurf_read_generic
479    tsurf(:,1) = tsurf_read_generic
480    tsoil(:,:,1) = tsoil_read_generic
481    emis(:,1) = emis_read_generic
[3028]482    watercap(:,1) = 0.
483    watercaptag(:) = .false.
484    albedo(:,1,1) = albedo_read_generic(:,1)
485    albedo(:,2,1) = albedo_read_generic(:,2)
[3149]486    inertiesoil(:,:,1) = inertiedat
[2842]487
[3028]488    if (nslope == 1) then
489        def_slope(1) = 0
490        def_slope(2) = 0
491        def_slope_mean = 0
492        subslope_dist(:,1) = 1.
493    endif
[2842]494
[3070]495    ! Remove unphysical values of surface tracer
[3149]496    qsurf(:,:,1) = qsurf_read_generic
[3070]497    where (qsurf < 0.) qsurf = 0.
[3591]498
499    deallocate(tsurf_read_generic,qsurf_read_generic,qsoil_read_generic,emis_read_generic)
[2842]500#endif
501
[3603]502do nnq = 1,nqtot  ! Why not using ini_tracer?
[3028]503    if (noms(nnq) == "h2o_ice") igcm_h2o_ice = nnq
504    if (noms(nnq) == "h2o_vap") then
505        igcm_h2o_vap = nnq
[3143]506        mmol(igcm_h2o_vap) = 18.
[3028]507    endif
508    if (noms(nnq) == "co2") igcm_co2 = nnq
[3065]509enddo
[3028]510
[2835]511!------------------------
[3028]512! I   Initialization
[2835]513!    I_c Subslope parametrisation
514!------------------------
[3028]515! Define some slope statistics
516iflat = 1
517do islope = 2,nslope
518    if (abs(def_slope_mean(islope)) < abs(def_slope_mean(iflat))) iflat = islope
519enddo
520write(*,*) 'Flat slope for islope = ',iflat
[3603]521write(*,*) 'Corresponding criterium = ',def_slope_mean(iflat)
[2794]522
523!------------------------
[3028]524! I   Initialization
[3161]525!    I_d Read the PCM data and convert them to the physical grid
[3028]526!------------------------
[3096]527! First we read the evolution of water and co2 ice (and the mass mixing ratio) over the first year of the PCM run, saving only the minimum value
[3571]528call get_timelen_PCM("data_PCM_Y1.nc",timelen)
[2794]529
[3591]530allocate(vmr_co2_PCM(ngrid,timelen),q_co2_PEM_phys(ngrid,timelen),q_h2o_PEM_phys(ngrid,timelen))
531allocate(ps_timeseries(ngrid,timelen),ps_avg(ngrid))
[3553]532allocate(min_co2_ice(ngrid,nslope,2),min_h2o_ice(ngrid,nslope,2))
[3591]533allocate(tsurf_avg_yr1(ngrid,nslope),tsurf_avg(ngrid,nslope))
[3599]534allocate(tsoil_avg(ngrid,nsoilmx,nslope),tsoil_timeseries(ngrid,nsoilmx,nslope,timelen))
[3591]535allocate(watersurf_density_avg(ngrid,nslope),watersoil_density_timeseries(ngrid,nsoilmx,nslope,timelen))
[2794]536
[3571]537call read_data_PCM("data_PCM_Y1.nc","data_PCM_Y2.nc",timelen,iim,jjm_value,ngrid,nslope,vmr_co2_PCM,ps_timeseries,ps_avg,tsurf_avg_yr1,tsurf_avg, &
538                   tsoil_avg,tsoil_timeseries,min_co2_ice,min_h2o_ice,q_co2_PEM_phys,q_h2o_PEM_phys,watersurf_density_avg,watersoil_density_timeseries)
[2985]539
[3571]540! Compute the deviation from the average
[3591]541allocate(ps_dev(ngrid),tsurf_dev(ngrid,nslope),tsoil_dev(ngrid,nsoilmx,nslope))
542ps_dev = ps_start0 - ps_avg
[3571]543tsurf_dev = tsurf - tsurf_avg
544tsoil_dev = tsoil - tsoil_avg(:,1:nsoilmx,:)
[2794]545
[2835]546!------------------------
[3028]547! I   Initialization
548!    I_e Initialization of the PEM variables and soil
[2835]549!------------------------
[3028]550call end_comsoil_h_PEM
551call ini_comsoil_h_PEM(ngrid,nslope)
552call end_adsorption_h_PEM
553call ini_adsorption_h_PEM(ngrid,nslope,nsoilmx_PEM)
[3493]554call end_ice_table
555call ini_ice_table(ngrid,nslope,nsoilmx_PEM)
[2794]556
[3599]557allocate(tsoil_PEM_timeseries(ngrid,nsoilmx_PEM,nslope,timelen),watersoil_density_PEM_avg(ngrid,nsoilmx_PEM,nslope),watersoil_density_PEM_timeseries(ngrid,nsoilmx_PEM,nslope,timelen))
[3028]558if (soil_pem) then
559    call soil_settings_PEM(ngrid,nslope,nsoilmx_PEM,nsoilmx,inertiesoil,TI_PEM)
[3367]560    tsoil_PEM(:,1:nsoilmx,:) = tsoil_avg
[3149]561    watersoil_density_PEM_timeseries(:,1:nsoilmx,:,:) = watersoil_density_timeseries
[3571]562    tsoil_PEM_timeseries(:,1:nsoilmx,:,:) = tsoil_timeseries
[3070]563    do l = nsoilmx + 1,nsoilmx_PEM
[3367]564        tsoil_PEM(:,l,:) = tsoil_avg(:,nsoilmx,:)
[3070]565        watersoil_density_PEM_timeseries(:,l,:,:) = watersoil_density_timeseries(:,nsoilmx,:,:)
[3571]566        tsoil_PEM_timeseries(:,l,:,:) = tsoil_timeseries(:,nsoilmx,:,:)
[3028]567    enddo
[3367]568    watersoil_density_PEM_avg = sum(watersoil_density_PEM_timeseries,4)/timelen
[3028]569endif !soil_pem
[3571]570deallocate(tsoil_avg,watersoil_density_timeseries,tsoil_timeseries)
[2794]571
[2779]572!------------------------
[3028]573! I   Initialization
[3143]574!    I_f Compute tendencies
[3028]575!------------------------
[3591]576allocate(d_co2ice(ngrid,nslope),d_h2oice(ngrid,nslope),d_co2ice_ini(ngrid,nslope))
[3498]577call compute_tend(ngrid,nslope,min_co2_ice,d_co2ice)
578call compute_tend(ngrid,nslope,min_h2o_ice,d_h2oice)
579d_co2ice_ini = d_co2ice
[3365]580deallocate(min_co2_ice,min_h2o_ice)
[2895]581
[2835]582!------------------------
[3028]583! I   Initialization
[3571]584!    I_g Compute global surface pressure
[3028]585!------------------------
[3571]586total_surface = sum(cell_area)
587ps_avg_global_ini = sum(cell_area*ps_avg)/total_surface
588ps_avg_global_new = ps_avg_global_ini
[3673]589write(*,*) "Total surface of the planet     =", total_surface
590write(*,*) "Initial global average pressure =", ps_avg_global_ini
[3384]591
592!------------------------
593! I   Initialization
594!    I_h Read the "startpem.nc"
595!------------------------
[3603]596write(*,*) '> Reading "startpem.nc"'
[3778]597allocate(co2_ice(ngrid,nslope),h2o_ice(ngrid,nslope),layerings_map(ngrid,nslope))
[3591]598allocate(delta_h2o_adsorbed(ngrid),delta_co2_adsorbed(ngrid),delta_h2o_icetablesublim(ngrid))
599delta_h2o_icetablesublim = 0.
[3602]600call pemetat0("startpem.nc",ngrid,nsoilmx,nsoilmx_PEM,nslope,timelen,timestep,TI_PEM,tsoil_PEM,icetable_depth,icetable_thickness,ice_porefilling, &
601              tsurf_avg_yr1,tsurf_avg,q_co2_PEM_phys,q_h2o_PEM_phys,ps_timeseries,ps_avg_global_ini,d_h2oice,d_co2ice,co2_ice,h2o_ice,            &
[3778]602              watersurf_density_avg,watersoil_density_PEM_avg,co2_adsorbed_phys,delta_co2_adsorbed,h2o_adsorbed_phys,delta_h2o_adsorbed,layerings_map)
[3571]603deallocate(tsurf_avg_yr1)
[2779]604
[3384]605! We save the places where h2o ice is sublimating
606! We compute the surface of h2o ice sublimating
[3591]607allocate(is_co2ice_sublim_ini(ngrid,nslope),is_h2oice_sublim_ini(ngrid,nslope),is_co2ice_ini(ngrid,nslope),co2ice_disappeared(ngrid,nslope))
[3571]608co2ice_sublim_surf_ini = 0.
[3384]609h2oice_ini_surf = 0.
[3571]610is_co2ice_sublim_ini = .false.
611is_h2oice_sublim_ini = .false.
[3384]612is_co2ice_ini = .false.
613co2ice_disappeared = .false.
[3779]614totmass_co2ice_ini = 0.
615totmass_atmco2_ini = 0.
[3770]616if (layering_algo) then
617    do ig = 1,ngrid
618        do islope = 1,nslope
[3809]619            if (is_co2ice_str(layerings_map(ig,islope)%top)) then
[3938]620                co2_ice(ig,islope) = layerings_map(ig,islope)%top%h_co2ice*rho_co2ice
[3809]621            else if (is_h2oice_str(layerings_map(ig,islope)%top)) then
[3938]622                h2o_ice(ig,islope) = layerings_map(ig,islope)%top%h_h2oice*rho_h2oice
[3809]623            else
624                call subsurface_ice_layering(layerings_map(ig,islope),h2o_ice_depth(ig,islope),h2o_ice(ig,islope),co2_ice(ig,islope))
625            endif
[3770]626        enddo
627    enddo
[3782]628    ! We put the sublimating tendency coming from subsurface ice into the overall tendency
[3842]629    where (h2o_ice_depth > 0. .and. zdqsdif_ssi_tot < -1.e-10) d_h2oice = zdqsdif_ssi_tot
[3770]630endif
[3384]631do i = 1,ngrid
[3779]632    totmass_atmco2_ini = totmass_atmco2_ini + cell_area(i)*ps_avg(i)/g
[3384]633    do islope = 1,nslope
[3779]634        totmass_co2ice_ini = totmass_co2ice_ini + co2_ice(i,islope)*cell_area(i)*subslope_dist(i,islope)/cos(pi*def_slope_mean(islope)/180.)
[3384]635        if (co2_ice(i,islope) > 0.) is_co2ice_ini(i,islope) = .true.
[3498]636        if (d_co2ice(i,islope) < 0. .and. co2_ice(i,islope) > 0.) then
[3571]637            is_co2ice_sublim_ini(i,islope) = .true.
638            co2ice_sublim_surf_ini = co2ice_sublim_surf_ini + cell_area(i)*subslope_dist(i,islope)
[3384]639        endif
[3782]640        if (d_h2oice(i,islope) < 0.) then
641            if (h2o_ice(i,islope) > 0.) then
642                is_h2oice_sublim_ini(i,islope) = .true.
643                h2oice_ini_surf = h2oice_ini_surf + cell_area(i)*subslope_dist(i,islope)
644            else if (h2o_ice_depth(i,islope) > 0.) then
645                is_h2oice_sublim_ini(i,islope) = .true.
646            endif
[3384]647        endif
648    enddo
649enddo
[3603]650write(*,*) "Total initial surface of CO2 ice sublimating =", co2ice_sublim_surf_ini
651write(*,*) "Total initial surface of H2O ice sublimating =", h2oice_ini_surf
[3384]652
[3779]653totmass_adsco2_ini = 0.
654totmass_adsh2o = 0.
[3028]655if (adsorption_pem) then
656    do ig = 1,ngrid
[3070]657        do islope = 1,nslope
[3028]658            do l = 1,nsoilmx_PEM - 1
[3571]659                if (l == 1) then
[3779]660                   totmass_adsco2_ini = totmass_adsco2_ini + co2_adsorbed_phys(ig,l,islope)*(layer_PEM(l))* &
[3028]661                                       subslope_dist(ig,islope)/cos(pi*def_slope_mean(islope)/180.)*cell_area(ig)
[3779]662                   totmass_adsh2o = totmass_adsh2o + h2o_adsorbed_phys(ig,l,islope)*(layer_PEM(l))* &
[3028]663                                       subslope_dist(ig,islope)/cos(pi*def_slope_mean(islope)/180.)*cell_area(ig)
[3264]664                else
[3779]665                   totmass_adsco2_ini = totmass_adsco2_ini + co2_adsorbed_phys(ig,l,islope)*(layer_PEM(l) - layer_PEM(l-1))* &
[3264]666                                       subslope_dist(ig,islope)/cos(pi*def_slope_mean(islope)/180.)*cell_area(ig)
[3779]667                   totmass_adsh2o = totmass_adsh2o + h2o_adsorbed_phys(ig,l,islope)*(layer_PEM(l) - layer_PEM(l-1))* &
[3264]668                                       subslope_dist(ig,islope)/cos(pi*def_slope_mean(islope)/180.)*cell_area(ig)
669                endif
[3028]670            enddo
[2961]671        enddo
[3028]672    enddo
[3779]673    totmass_adsco2 = totmass_adsco2_ini
674    write(*,*) "Tot mass of CO2 in the regolith =", totmass_adsco2
675    write(*,*) "Tot mass of H2O in the regolith =", totmass_adsh2o
[3028]676endif ! adsorption
[2794]677
[2835]678!------------------------
[3028]679! I   Initialization
[2835]680!    I_i Compute orbit criterion
[3028]681!------------------------
[2842]682#ifndef CPP_STD
[3050]683    call iniorbit(aphelie,periheli,year_day,peri_day,obliquit)
[2842]684#else
[3050]685    call iniorbit(apoastr,periastr,year_day,peri_day,obliquit)
[2842]686#endif
[2794]687
[3498]688n_myear_leg = Max_iter_pem
689if (evol_orbit_pem) call orbit_param_criterion(i_myear,n_myear_leg)
[3403]690
[3028]691!-------------------------- END INITIALIZATION -------------------------
[2794]692
[3028]693!-------------------------------- RUN ----------------------------------
[2794]694!------------------------
695! II  Run
[3065]696!    II_a Update pressure, ice and tracers
[2794]697!------------------------
[3616]698write(*,*)
[3673]699write(*,*) '********* PEM cycle *********'
[3498]700i_myear_leg = 0
[3149]701stopPEM = 0
[3319]702if (layering_algo) then
[3784]703    allocate(h2o_ice_depth_old(ngrid,nslope),new_str(ngrid,nslope),new_lag(ngrid,nslope),current(ngrid,nslope))
[3319]704    new_str = .true.
[3770]705    new_lag = .true.
[3319]706    do islope = 1,nslope
707        do ig = 1,ngrid
[3778]708            current(ig,islope)%p => layerings_map(ig,islope)%top
[3319]709        enddo
[3297]710    enddo
[3319]711endif
[2794]712
[3498]713do while (i_myear_leg < n_myear_leg .and. i_myear < n_myear)
[2835]714! II.a.1. Compute updated global pressure
[3673]715    write(*,'(a,f10.2)') ' **** Iteration of the PEM leg (Martian years): ', i_myear_leg + 1
[3603]716    write(*,*) "> Updating the surface pressure"
[3571]717    ps_avg_global_old = ps_avg_global_new
[3028]718    do i = 1,ngrid
719        do islope = 1,nslope
[3611]720            ps_avg_global_new = ps_avg_global_new - CO2cond_ps*g*cell_area(i)*d_co2ice(i,islope)*subslope_dist(i,islope)/cos(pi*def_slope_mean(islope)/180.)/total_surface
[3028]721        enddo
722    enddo
723    if (adsorption_pem) then
724        do i = 1,ngrid
[3611]725            ps_avg_global_new = ps_avg_global_new - g*cell_area(i)*delta_co2_adsorbed(i)/total_surface
[3050]726        enddo
[3028]727    endif
[3680]728    ps_avg = ps_avg*ps_avg_global_new/ps_avg_global_old
[3603]729    write(*,*) 'Global average pressure old time step:',ps_avg_global_old
730    write(*,*) 'Global average pressure new time step:',ps_avg_global_new
[2835]731
[3571]732! II.a.2. Pressure timeseries (the values are deleted when unused because of big memory consumption)
[3603]733    write(*,*) "> Updating the surface pressure timeseries for the new pressure"
[3571]734    allocate(zplev_timeseries_old(ngrid,nlayer + 1,timelen))
[3028]735    do l = 1,nlayer + 1
736        do ig = 1,ngrid
[3571]737            zplev_timeseries_old(ig,l,:) = ap(l) + bp(l)*ps_timeseries(ig,:)
[3028]738        enddo
739    enddo
[3673]740    ps_timeseries(:,:) = ps_timeseries(:,:)*ps_avg_global_new/ps_avg_global_old
[3603]741    write(*,*) "> Updating the pressure levels timeseries for the new pressure"
[3571]742    allocate(zplev_timeseries_new(ngrid,nlayer + 1,timelen))
[3028]743    do l = 1,nlayer + 1
744        do ig = 1,ngrid
[3571]745            zplev_timeseries_new(ig,l,:) = ap(l) + bp(l)*ps_timeseries(ig,:)
[3028]746        enddo
747    enddo
[2779]748
[3571]749! II.a.3. Tracers timeseries
[3603]750    write(*,*) "> Updating the tracer VMR timeseries for the new pressure"
[3591]751    allocate(vmr_co2_PEM_phys(ngrid,timelen))
[3028]752    l = 1
753    do ig = 1,ngrid
754        do t = 1,timelen
[3571]755            ! H2O
756            q_h2o_PEM_phys(ig,t) = q_h2o_PEM_phys(ig,t)*(zplev_timeseries_old(ig,l,t) - zplev_timeseries_old(ig,l + 1,t))/ &
757                                   (zplev_timeseries_new(ig,l,t) - zplev_timeseries_new(ig,l + 1,t))
[3143]758            if (q_h2o_PEM_phys(ig,t) < 0) then
759                q_h2o_PEM_phys(ig,t) = 1.e-30
760            else if (q_h2o_PEM_phys(ig,t) > 1) then
761                q_h2o_PEM_phys(ig,t) = 1.
762            endif
[3571]763            ! CO2
764            q_co2_PEM_phys(ig,t) = q_co2_PEM_phys(ig,t)*(zplev_timeseries_old(ig,l,t) - zplev_timeseries_old(ig,l + 1,t))/ &
765                                   (zplev_timeseries_new(ig,l,t) - zplev_timeseries_new(ig,l + 1,t))                       &
766                                + ((zplev_timeseries_new(ig,l,t) - zplev_timeseries_new(ig,l + 1,t))                       &
767                                -  (zplev_timeseries_old(ig,l,t) - zplev_timeseries_old(ig,l + 1,t)))/                     &
768                                   (zplev_timeseries_new(ig,l,t) - zplev_timeseries_new(ig,l + 1,t))
[3028]769            if (q_co2_PEM_phys(ig,t) < 0) then
770                q_co2_PEM_phys(ig,t) = 1.e-30
[3143]771            else if (q_co2_PEM_phys(ig,t) > 1) then
[3028]772                q_co2_PEM_phys(ig,t) = 1.
773            endif
[3571]774            mmean = 1./(A*q_co2_PEM_phys(ig,t) + B)
[3149]775            vmr_co2_PEM_phys(ig,t) = q_co2_PEM_phys(ig,t)*mmean/m_co2
[3028]776        enddo
777    enddo
[3571]778    deallocate(zplev_timeseries_new,zplev_timeseries_old)
[2794]779
[3028]780!------------------------
[2835]781! II  Run
[3149]782!    II_b Evolution of ice
[3028]783!------------------------
[3591]784    allocate(zshift_surf(ngrid,nslope),zlag(ngrid,nslope))
[3319]785    if (layering_algo) then
[3784]786        h2o_ice_depth_old = h2o_ice_depth
[3938]787
788        allocate(d_h2oice_new(ngrid,nslope))
789        call stopping_crit_h2o(ngrid,nslope,cell_area,delta_h2o_adsorbed,delta_h2o_icetablesublim,h2o_ice,d_h2oice,S_atm_2_h2o,S_h2o_2_atm,S_atm_2_h2oice,S_h2oice_2_atm,stopPEM)
790        call balance_h2oice_reservoirs(ngrid,nslope,S_atm_2_h2o,S_h2o_2_atm,S_atm_2_h2oice,S_h2oice_2_atm,h2o_ice,d_h2oice,d_h2oice_new)
791   
[3319]792        do islope = 1,nslope
793            do ig = 1,ngrid
[3938]794                call make_layering(layerings_map(ig,islope),d_co2ice(ig,islope),d_h2oice_new(ig,islope),new_str(ig,islope),zshift_surf(ig,islope),new_lag(ig,islope),zlag(ig,islope),current(ig,islope)%p)
[3842]795                !call print_layering(layerings_map(ig,islope))
[3784]796                co2_ice(ig,islope) = 0.
797                h2o_ice(ig,islope) = 0.
798                h2o_ice_depth(ig,islope) = 0.
[3778]799                if (is_co2ice_str(layerings_map(ig,islope)%top)) then
[3938]800                    co2_ice(ig,islope) = layerings_map(ig,islope)%top%h_co2ice*rho_co2ice
[3778]801                else if (is_h2oice_str(layerings_map(ig,islope)%top)) then
[3938]802                    h2o_ice(ig,islope) = layerings_map(ig,islope)%top%h_h2oice*rho_h2oice
[3778]803                else
[3809]804                    call subsurface_ice_layering(layerings_map(ig,islope),h2o_ice_depth(ig,islope),h2o_ice(ig,islope),co2_ice(ig,islope))
[3778]805                endif
[3319]806            enddo
[3297]807        enddo
[3938]808        deallocate(d_h2oice_new)
[3553]809    else
810        zlag = 0.
[3770]811        call evol_h2o_ice(ngrid,nslope,cell_area,delta_h2o_adsorbed,delta_h2o_icetablesublim,h2o_ice,d_h2oice,zshift_surf,stopPEM)
812        call evol_co2_ice(ngrid,nslope,co2_ice,d_co2ice,zshift_surf)
[3319]813    endif
[2794]814
815!------------------------
816! II  Run
[3149]817!    II_c Flow of glaciers
[2794]818!------------------------
[3591]819    allocate(flag_co2flow(ngrid,nslope),flag_h2oflow(ngrid,nslope))
[3809]820    if (co2ice_flow .and. nslope > 1) then
821        call flow_co2glaciers(timelen,ngrid,nslope,iflat,subslope_dist,def_slope_mean,vmr_co2_PEM_phys, &
822                              ps_timeseries,ps_avg_global_old,ps_avg_global_new,co2_ice,flag_co2flow)
823        if (layering_algo) then
824            do islope = 1,nslope
825                do ig = 1,ngrid
[3938]826                    layerings_map(ig,islope)%top%h_co2ice = co2_ice(ig,islope)/rho_co2ice
[3809]827                enddo
[3770]828            enddo
[3809]829        endif
[3770]830    endif
[3809]831    if (h2oice_flow .and. nslope > 1) then
832        call flow_h2oglaciers(ngrid,nslope,iflat,subslope_dist,def_slope_mean,tsurf_avg,h2o_ice,flag_h2oflow)
833        if (layering_algo) then
834            do islope = 1,nslope
835                do ig = 1,ngrid
[3938]836                    layerings_map(ig,islope)%top%h_h2oice = h2o_ice(ig,islope)/rho_h2oice
[3809]837                enddo
838            enddo
839        endif
840    endif
[3065]841
[2794]842!------------------------
843! II  Run
[2835]844!    II_d Update surface and soil temperatures
[2794]845!------------------------
[2835]846! II_d.1 Update Tsurf
[3933]847    call update_tsurf_nearest_baresoil(ngrid,nslope,iim,jjm_value,tsurf_avg,co2_ice,is_co2ice_ini,co2ice_disappeared)
[2794]848
[3028]849    if (soil_pem) then
[3553]850! II_d.2 Shifting soil temperature to surface
[3770]851        call shift_tsoil2surf(ngrid,nsoilmx_PEM,nslope,zshift_surf,zlag,tsurf_avg,tsoil_PEM)
[2794]852
[3553]853! II_d.3 Update soil temperature
[3603]854        write(*,*)"> Updating soil temperature profile"
[3610]855        allocate(tsoil_avg_old(ngrid,nsoilmx_PEM),tsoil_PEM_timeseries_old(ngrid,nsoilmx_PEM,nslope,timelen))
856        tsoil_PEM_timeseries_old = tsoil_PEM_timeseries
[3426]857        do islope = 1,nslope
[3571]858            tsoil_avg_old = tsoil_PEM(:,:,islope)
[3426]859            call compute_tsoil_pem(ngrid,nsoilmx_PEM,.true.,TI_PEM(:,:,islope),timestep,tsurf_avg(:,islope),tsoil_PEM(:,:,islope))
860            call compute_tsoil_pem(ngrid,nsoilmx_PEM,.false.,TI_PEM(:,:,islope),timestep,tsurf_avg(:,islope),tsoil_PEM(:,:,islope))
[2794]861
[3028]862            do t = 1,timelen
863                do ig = 1,ngrid
864                    do isoil = 1,nsoilmx_PEM
[3609]865                        ! Update of soil temperature timeseries which is needed to compute the water soil density timeseries
866                        tsoil_PEM_timeseries(ig,isoil,islope,t) = tsoil_PEM_timeseries(ig,isoil,islope,t)*tsoil_PEM(ig,isoil,islope)/tsoil_avg_old(ig,isoil)
867                        ! Update of watersoil density
[3571]868                        watersoil_density_PEM_timeseries(ig,isoil,islope,t) = exp(beta_clap_h2o/tsoil_PEM_timeseries(ig,isoil,islope,t) + alpha_clap_h2o)/tsoil_PEM_timeseries(ig,isoil,islope,t)*mmol(igcm_h2o_vap)/(mugaz*r)
[3426]869                        if (isnan(tsoil_PEM(ig,isoil,islope))) call abort_pem("PEM - Update Tsoil","NaN detected in tsoil_PEM",1)
[3028]870                    enddo
871                enddo
872            enddo
873        enddo
[3367]874        watersoil_density_PEM_avg = sum(watersoil_density_PEM_timeseries,4)/timelen
[3571]875        deallocate(tsoil_avg_old)
[2888]876
[3553]877! II_d.4 Update the ice table
[3610]878        allocate(icetable_thickness_old(ngrid,nslope),ice_porefilling_old(ngrid,nsoilmx_PEM,nslope),icetable_depth_old(ngrid,nslope))
[3170]879        if (icetable_equilibrium) then
[3603]880            write(*,*) "> Updating ice table (equilibrium method)"
[3525]881            icetable_thickness_old = icetable_thickness
[3493]882            call computeice_table_equilibrium(ngrid,nslope,nsoilmx_PEM,watercaptag,watersurf_density_avg,watersoil_density_PEM_avg,TI_PEM(:,1,:),icetable_depth,icetable_thickness)
[3525]883            call compute_massh2o_exchange_ssi(ngrid,nslope,nsoilmx_PEM,icetable_thickness_old,ice_porefilling_old,tsurf_avg,tsoil_PEM,delta_h2o_icetablesublim) ! Mass of H2O exchange between the ssi and the atmosphere
[3490]884        else if (icetable_dynamic) then
[3603]885            write(*,*) "> Updating ice table (dynamic method)"
[3525]886            ice_porefilling_old = ice_porefilling
[3610]887            icetable_depth_old = icetable_depth
[3493]888            allocate(porefill(nsoilmx_PEM))
889            do ig = 1,ngrid
890                do islope = 1,nslope
[3571]891                    call dyn_ss_ice_m(icetable_depth(ig,islope),tsurf_avg(ig,islope),tsoil_PEM(ig,:,islope),nsoilmx_PEM,TI_PEM(ig,1,nslope),ps_avg(ig),(/sum(q_h2o_PEM_phys(ig,:))/size(q_h2o_PEM_phys,2)/),ice_porefilling(ig,:,islope),porefill,ssi_depth)
[3493]892                    icetable_depth(ig,islope) = ssi_depth
893                    ice_porefilling(ig,:,islope) = porefill
894                enddo
895            enddo
896            deallocate(porefill)
[3673]897            call compute_massh2o_exchange_ssi(ngrid,nslope,nsoilmx_PEM,icetable_thickness_old,ice_porefilling_old,tsurf_avg,tsoil_PEM,delta_h2o_icetablesublim) ! Mass of H2O exchange between the ssi and the atmosphere
[3486]898        endif
[3591]899        deallocate(icetable_thickness_old,ice_porefilling_old)
[3486]900
[3553]901! II_d.5 Update the soil thermal properties
[3571]902        call update_soil_thermalproperties(ngrid,nslope,nsoilmx_PEM,d_h2oice,h2o_ice,ps_avg_global_new,icetable_depth,icetable_thickness,ice_porefilling,icetable_equilibrium,icetable_dynamic,TI_PEM)
[2794]903
[3553]904! II_d.6 Update the mass of the regolith adsorbed
[3779]905        totmass_adsco2 = 0.
906        totmass_adsh2o = 0.
[3028]907        if (adsorption_pem) then
[3571]908            call regolith_adsorption(ngrid,nslope,nsoilmx_PEM,timelen,d_h2oice,d_co2ice,h2o_ice,co2_ice, &
909                                     tsoil_PEM,TI_PEM,ps_timeseries,q_co2_PEM_phys,q_h2o_PEM_phys,       &
[3553]910                                     h2o_adsorbed_phys,delta_h2o_adsorbed,co2_adsorbed_phys,delta_co2_adsorbed)
[3028]911            do ig = 1,ngrid
[3493]912                do islope = 1,nslope
[3532]913                    do l = 1,nsoilmx_PEM
[3493]914                        if (l == 1) then
[3779]915                            totmass_adsco2 = totmass_adsco2 + co2_adsorbed_phys(ig,l,islope)*(layer_PEM(l))* &
[3264]916                                       subslope_dist(ig,islope)/cos(pi*def_slope_mean(islope)/180.)*cell_area(ig)
[3779]917                            totmass_adsh2o = totmass_adsh2o + h2o_adsorbed_phys(ig,l,islope)*(layer_PEM(l))* &
[3264]918                                       subslope_dist(ig,islope)/cos(pi*def_slope_mean(islope)/180.)*cell_area(ig)
919                        else
[3779]920                            totmass_adsco2 = totmass_adsco2 + co2_adsorbed_phys(ig,l,islope)*(layer_PEM(l) - layer_PEM(l - 1))* &
[3264]921                                       subslope_dist(ig,islope)/cos(pi*def_slope_mean(islope)/180.)*cell_area(ig)
[3779]922                            totmass_adsh2o = totmass_adsh2o + h2o_adsorbed_phys(ig,l,islope)*(layer_PEM(l) - layer_PEM(l - 1))* &
[3264]923                                       subslope_dist(ig,islope)/cos(pi*def_slope_mean(islope)/180.)*cell_area(ig)
924                        endif
[3028]925                    enddo
926                enddo
927            enddo
[3779]928            write(*,*) "Total mass of CO2 in the regolith =", totmass_adsco2
929            write(*,*) "Total mass of H2O in the regolith =", totmass_adsh2o
[3028]930        endif
931    endif !soil_pem
[3666]932    deallocate(zshift_surf,zlag)
[3028]933
[2794]934!------------------------
935! II  Run
[3088]936!    II_e Outputs
[2794]937!------------------------
[3571]938    call writediagpem(ngrid,'ps_avg','Global average pressure','Pa',0,(/ps_avg_global_new/))
[3088]939    do islope = 1,nslope
940        write(str2(1:2),'(i2.2)') islope
[3181]941        call writediagpem(ngrid,'h2o_ice_slope'//str2,'H2O ice','kg.m-2',2,h2o_ice(:,islope))
942        call writediagpem(ngrid,'co2_ice_slope'//str2,'CO2 ice','kg.m-2',2,co2_ice(:,islope))
[3498]943        call writediagpem(ngrid,'d_h2oice_slope'//str2,'H2O ice tend','kg.m-2.year-1',2,d_h2oice(:,islope))
944        call writediagpem(ngrid,'d_co2ice_slope'//str2,'CO2 ice tend','kg.m-2.year-1',2,d_co2ice(:,islope))
[3591]945        call writediagpem(ngrid,'Flow_co2ice_slope'//str2,'CO2 ice flow','Boolean',2,real(flag_co2flow(:,islope)))
946        call writediagpem(ngrid,'Flow_h2oice_slope'//str2,'H2O ice flow','Boolean',2,real(flag_h2oflow(:,islope)))
[3571]947        call writediagpem(ngrid,'tsurf_slope'//str2,'tsurf','K',2,tsurf_avg(:,islope))
[3339]948        if (icetable_equilibrium) then
[3493]949            call writediagpem(ngrid,'ssi_depth_slope'//str2,'ice table depth','m',2,icetable_depth(:,islope))
[3537]950            call writediagpem(ngrid,'ssi_thick_slope'//str2,'ice table thickness','m',2,icetable_thickness(:,islope))
[3490]951        else if (icetable_dynamic) then
[3493]952            call writediagpem(ngrid,'ssi_depth_slope'//str2,'ice table depth','m',2,icetable_depth(:,islope))
[3486]953        endif
954
[3339]955        if (soil_pem) then
[3571]956            call writediagsoilpem(ngrid,'tsoil_PEM_slope'//str2,'tsoil','K',3,tsoil_PEM(:,:,islope))
957            call writediagsoilpem(ngrid,'inertiesoil_PEM_slope'//str2,'TI','K',3,TI_PEM(:,:,islope))
[3543]958            if (icetable_dynamic) call writediagsoilpem(ngrid,'ice_porefilling'//str2,'ice pore filling','-',3,ice_porefilling(:,:,islope))
[3171]959            if (adsorption_pem) then
[3553]960                call writediagsoilpem(ngrid,'co2_ads_slope'//str2,'co2_ads','K',3,co2_adsorbed_phys(:,:,islope))
961                call writediagsoilpem(ngrid,'h2o_ads_slope'//str2,'h2o_ads','K',3,h2o_adsorbed_phys(:,:,islope))
[3532]962            endif
[3171]963        endif
[3088]964    enddo
[3591]965    deallocate(flag_co2flow,flag_h2oflow)
[3088]966
[3779]967    ! Checking mass balance for CO2
[3851]968    if (abs(CO2cond_ps - 1.) < 1.e-10) then
969        totmass_co2ice = 0.
970        totmass_atmco2 = 0.
971        do ig = 1,ngrid
972            totmass_atmco2 = totmass_atmco2 + cell_area(ig)*ps_avg(ig)/g
973            do islope = 1,nslope
974                totmass_co2ice = totmass_co2ice + co2_ice(ig,islope)*cell_area(ig)*subslope_dist(ig,islope)/cos(pi*def_slope_mean(islope)/180.)
975            enddo
[3842]976        enddo
[3851]977        totmass_ini = max(totmass_atmco2_ini + totmass_co2ice_ini + totmass_adsco2_ini,1.e-10)
978        write(*,'(a,f8.3,a)') " > Relative total CO2 mass balance = ", 100.*(totmass_atmco2 + totmass_co2ice + totmass_adsco2 - totmass_atmco2_ini - totmass_co2ice_ini - totmass_adsco2_ini)/totmass_ini, ' %'
[3938]979        if (abs((totmass_atmco2 + totmass_co2ice + totmass_adsco2 - totmass_atmco2_ini - totmass_co2ice_ini - totmass_adsco2_ini)/totmass_ini) > 0.01) then
980            write(*,*) '  /!\ Warning: mass balance is not conserved!'
[3851]981            totmass_ini = max(totmass_atmco2_ini,1.e-10)
[3938]982            write(*,*) '       Atmospheric CO2 mass balance = ', 100.*(totmass_atmco2 - totmass_atmco2_ini)/totmass_ini, ' %'
[3851]983            totmass_ini = max(totmass_co2ice_ini,1.e-10)
[3938]984            write(*,*) '       CO2 ice mass balance         = ', 100.*(totmass_co2ice - totmass_co2ice_ini)/totmass_ini, ' %'
[3851]985            totmass_ini = max(totmass_adsco2_ini,1.e-10)
[3938]986            write(*,*) '       Adsorbed CO2 mass balance    = ', 100.*(totmass_adsco2 - totmass_adsco2_ini)/totmass_ini, ' %'
[3851]987        endif
[3779]988    endif
989
[3088]990!------------------------
991! II  Run
992!    II_f Update the tendencies
993!------------------------
[3571]994    call recomp_tend_co2(ngrid,nslope,timelen,d_co2ice,d_co2ice_ini,co2_ice,emis,vmr_co2_PCM,vmr_co2_PEM_phys,ps_timeseries,ps_avg_global_old,ps_avg_global_new)
[3707]995    write(*,*) "> Updating the H2O sub-surface ice tendency due to lag layer"
[3770]996    if (layering_algo) then
997        do ig = 1,ngrid
998            do islope = 1,nslope
[3789]999                if (is_h2oice_sublim_ini(ig,islope) .and. h2o_ice_depth(ig,islope) > 0.) call recomp_tend_h2o(h2o_ice_depth_old(ig,islope),h2o_ice_depth(ig,islope),tsurf_avg(ig,islope),tsoil_PEM_timeseries_old(ig,:,islope,:),tsoil_PEM_timeseries(ig,:,islope,:),d_h2oice(ig,islope))
[3770]1000            enddo
[3707]1001        enddo
[3770]1002!~     else
1003!~         do ig = 1,ngrid
1004!~             do islope = 1,nslope
1005!~                 call recomp_tend_h2o(icetable_depth_old(ig,islope),icetable_depth(ig,islope),tsurf_avg(ig,islope),tsoil_PEM_timeseries_old(ig,:,islope,:),tsoil_PEM_timeseries(ig,:,islope,:),d_h2oice(ig,islope))
1006!~             enddo
1007!~         enddo
1008    endif
[3666]1009    if (soil_pem) deallocate(icetable_depth_old,tsoil_PEM_timeseries_old)
1010    deallocate(vmr_co2_PEM_phys)
[2794]1011
[2835]1012!------------------------
1013! II  Run
[3088]1014!    II_g Checking the stopping criterion
[2835]1015!------------------------
[3603]1016    write(*,*) "> Checking the stopping criteria"
[3571]1017    call stopping_crit_h2o_ice(cell_area,h2oice_ini_surf,is_h2oice_sublim_ini,h2o_ice,stopPEM,ngrid)
1018    call stopping_crit_co2(cell_area,co2ice_sublim_surf_ini,is_co2ice_sublim_ini,co2_ice,stopPEM,ngrid,ps_avg_global_ini,ps_avg_global_new,nslope)
[3498]1019    i_myear_leg = i_myear_leg + dt
1020    i_myear = i_myear + dt
1021    if (stopPEM <= 0 .and. i_myear_leg >= n_myear_leg) stopPEM = 5
[3430]1022    if (stopPEM <= 0 .and. i_myear >= n_myear) stopPEM = 6
[3389]1023    call system_clock(c2)
[3430]1024    if (stopPEM <= 0 .and. timewall .and. real((c2 - c1)/cr) >= timelimit - antetime) stopPEM = 7
[3432]1025    if (stopPEM > 0) then
[3149]1026        select case (stopPEM)
1027            case(1)
[3673]1028                write(*,'(a,i0,a)') " **** STOPPING because surface of h2o ice sublimating is too low: ", stopPEM, ". See message above."
[3149]1029            case(2)
[3673]1030                write(*,'(a,i0,a)') " **** STOPPING because tendencies on h2o ice = 0: ", stopPEM, ". See message above."
[3149]1031            case(3)
[3673]1032                write(*,'(a,i0,a)') " **** STOPPING because surface of co2 ice sublimating is too low: ", stopPEM, ". See message above."
[3149]1033            case(4)
[3673]1034                write(*,'(a,i0,a)') " **** STOPPING because surface global pressure changed too much: ", stopPEM, ". See message above."
[3149]1035            case(5)
[3673]1036                write(*,'(a,i0)') " **** STOPPING because maximum number of iterations is reached (possibly due to orbital parameters): ", stopPEM
[3149]1037            case(6)
[3673]1038                write(*,'(a,i0)') " **** STOPPING because maximum number of Martian years to be simulated is reached: ", stopPEM
[3363]1039            case(7)
[3673]1040                write(*,'(a,i0)') " **** STOPPING because the time limit for the PEM job will be reached soon: ", stopPEM
[3149]1041            case default
[3673]1042                write(*,'(a,i0)') " **** STOPPING with unknown stopping criterion code: ", stopPEM
[3149]1043        end select
[2779]1044        exit
[3028]1045    else
[3673]1046        write(*,'(a,f10.2,a)') ' **** The chained simulation has run for ',i_myear,' Martian years.'
1047        write(*,*) '**** The PEM can continue!'
1048        write(*,*) '****'
[3028]1049    endif
[3149]1050enddo ! big time iteration loop of the pem
[3591]1051deallocate(vmr_co2_PCM,q_co2_PEM_phys,q_h2o_PEM_phys,delta_co2_adsorbed)
[3598]1052deallocate(watersoil_density_PEM_avg,watersurf_density_avg)
[3591]1053deallocate(ps_timeseries,tsoil_PEM_timeseries,watersoil_density_PEM_timeseries)
1054deallocate(co2ice_disappeared,delta_h2o_adsorbed,delta_h2o_icetablesublim)
1055deallocate(d_co2ice,d_co2ice_ini,d_h2oice)
1056deallocate(is_co2ice_ini,is_co2ice_sublim_ini,is_h2oice_sublim_ini)
[3784]1057if (layering_algo) deallocate(h2o_ice_depth_old,new_str,new_lag,current)
[3028]1058!------------------------------ END RUN --------------------------------
[2779]1059
[3028]1060!------------------------------- OUTPUT --------------------------------
[2794]1061!------------------------
1062! III Output
[3770]1063!    III_a Update surface values for the PCM start files
[2794]1064!------------------------
[3616]1065write(*,*)
[3673]1066write(*,*) '********* PEM finalization *********'
[3778]1067! III_a.1 Ice update for start file
[3603]1068write(*,*) '> Reconstructing perennial ice and frost for the PCM'
[3149]1069watercap = 0.
[3770]1070perennial_co2ice = co2_ice
[3028]1071do ig = 1,ngrid
[3159]1072    ! H2O ice metamorphism
[3705]1073    !if (metam_h2oice .and. sum(qsurf(ig,igcm_h2o_ice,:)*subslope_dist(ig,:)/cos(pi*def_slope_mean(:)/180.)) > metam_h2oice_threshold) then
1074    !    h2o_ice(ig,:) = h2o_ice(ig,:) + qsurf(ig,igcm_h2o_ice,:) - metam_h2oice_threshold
1075    !    qsurf(ig,igcm_h2o_ice,:) = metam_h2oice_threshold
1076    !endif
[3159]1077
1078    ! Is H2O ice still considered as an infinite reservoir for the PCM?
[3770]1079    if (sum(h2o_ice(ig,:)*subslope_dist(ig,:)/cos(pi*def_slope_mean(:)/180.)) > inf_h2oice_threshold) then
1080        ! There is enough ice to be considered as an infinite reservoir
[3149]1081        watercaptag(ig) = .true.
1082    else
[3770]1083        ! Too little ice to be considered as an infinite reservoir so ice is transferred to the frost
[3149]1084        watercaptag(ig) = .false.
[3770]1085        qsurf(ig,igcm_h2o_ice,:) = qsurf(ig,igcm_h2o_ice,:) + h2o_ice(ig,:)
[3149]1086        h2o_ice(ig,:) = 0.
[3028]1087    endif
[3159]1088
1089    ! CO2 ice metamorphism
[3705]1090    !if (metam_co2ice .and. sum(qsurf(ig,igcm_co2,:)*subslope_dist(ig,:)/cos(pi*def_slope_mean(:)/180.)) > metam_co2ice_threshold) then
1091    !    perennial_co2ice(ig,:) = perennial_co2ice(ig,:) + qsurf(ig,igcm_co2,:) - metam_co2ice_threshold
1092    !    qsurf(ig,igcm_co2,:) = metam_co2ice_threshold
1093    !endif
[3028]1094enddo
[2888]1095
[3770]1096! III.a.3. Tsurf update for start file
[3603]1097write(*,*) '> Reconstructing the surface temperature for the PCM'
[3571]1098tsurf = tsurf_avg + tsurf_dev
[3591]1099deallocate(tsurf_dev)
[3571]1100
[3770]1101! III_a.4 Tsoil update for start file
[3028]1102if (soil_pem) then
[3603]1103    write(*,*) '> Reconstructing the soil temperature profile for the PCM'
[3532]1104    inertiesoil = TI_PEM(:,:nsoilmx,:)
[3571]1105    ! Tsurf has evolved and so the soil temperature profile needs to be adapted to match this new value
1106    do isoil = 1,nsoilmx
1107        tsoil_dev(:,isoil,:) = tsoil_dev(:,isoil,:)*(tsurf_avg(:,:) - tsoil_PEM(:,1,:))/tsoil_dev(:,1,:)
1108    enddo
1109    tsoil = tsoil_PEM(:,1:nsoilmx,:) + tsoil_dev
[3172]1110#ifndef CPP_STD
[3181]1111    flux_geo = fluxgeo
[3172]1112#endif
[3028]1113endif
[3591]1114deallocate(tsurf_avg,tsoil_dev)
[2779]1115
[3770]1116! III_a.5 Pressure update for start file
[3603]1117write(*,*) '> Reconstructing the pressure for the PCM'
[3591]1118allocate(ps_start(ngrid))
[3680]1119! The pressure deviation is rescaled as well to avoid disproportionate oscillations in case of huge average pressure drop
1120ps_start = ps_avg + ps_dev*ps_avg_global_new/ps_avg_global_ini
[3591]1121deallocate(ps_avg,ps_dev)
[2794]1122
[3770]1123! III_a.6 Tracers update for start file
[3603]1124write(*,*) '> Reconstructing the tracer VMR for the PCM'
[3571]1125allocate(zplev_start0(ngrid,nlayer + 1),zplev_new(ngrid,nlayer + 1))
[3028]1126do l = 1,nlayer + 1
[3571]1127    zplev_start0(:,l) = ap(l) + bp(l)*ps_start0
1128    zplev_new(:,l) = ap(l) + bp(l)*ps_start
[3028]1129enddo
[2835]1130
[3028]1131do nnq = 1,nqtot
1132    if (noms(nnq) /= "co2") then
1133        do l = 1,llm - 1
1134            do ig = 1,ngrid
[3571]1135                q(ig,l,nnq) = q(ig,l,nnq)*(zplev_start0(ig,l) - zplev_start0(ig,l + 1))/(zplev_new(ig,l) - zplev_new(ig,l + 1))
[3028]1136            enddo
1137            q(:,llm,nnq) = q(:,llm - 1,nnq)
1138        enddo
1139    else
1140        do l = 1,llm - 1
1141            do ig = 1,ngrid
[3571]1142                q(ig,l,nnq) = q(ig,l,nnq)*(zplev_start0(ig,l) - zplev_start0(ig,l + 1))/(zplev_new(ig,l) - zplev_new(ig,l + 1)) &
1143                              + ((zplev_new(ig,l) - zplev_new(ig,l + 1)) - (zplev_start0(ig,l) - zplev_start0(ig,l + 1)))/(zplev_new(ig,l) - zplev_new(ig,l + 1))
[3028]1144            enddo
1145            q(:,llm,nnq) = q(:,llm - 1,nnq)
1146        enddo
1147    endif
1148enddo
[3591]1149deallocate(zplev_start0)
[2835]1150
[3571]1151! Conserving the tracers mass for start file
[3028]1152do nnq = 1,nqtot
1153    do ig = 1,ngrid
1154        do l = 1,llm - 1
1155            if (q(ig,l,nnq) > 1 .and. (noms(nnq) /= "dust_number") .and. (noms(nnq) /= "ccn_number") .and. (noms(nnq) /= "stormdust_number") .and. (noms(nnq) /= "topdust_number")) then
[3065]1156                extra_mass = (q(ig,l,nnq) - 1)*(zplev_new(ig,l) - zplev_new(ig,l + 1))
1157                q(ig,l,nnq) = 1.
1158                q(ig,l + 1,nnq) = q(ig,l + 1,nnq) + extra_mass*(zplev_new(ig,l + 1) - zplev_new(ig,l + 2))
[3028]1159                write(*,*) 'extra ',noms(nnq),extra_mass, noms(nnq) /= "dust_number",noms(nnq) /= "ccn_number"
[3770]1160            endif
[3028]1161            if (q(ig,l,nnq) < 0) q(ig,l,nnq) = 1.e-30
1162        enddo
1163    enddo
1164enddo
[3591]1165deallocate(zplev_new)
[2779]1166
[3770]1167! III_a.7 Albedo update for start file
[3603]1168write(*,*) '> Reconstructing the albedo for the PCM'
[3584]1169do ig = 1,ngrid
1170    if (latitude(ig) < 0.) then
1171        icap = 2 ! Southern hemisphere
1172    else
1173        icap = 1 ! Northern hemisphere
1174    endif
[3603]1175    do islope = 1,nslope
[3584]1176        ! Bare ground
1177        albedo(ig,:,islope) = albedodat(ig)
1178        emis(ig,islope) = emissiv
1179
1180        ! CO2 ice/frost is treated after H20 ice/frost because it is considered dominant
1181        ! H2O ice
1182        if (h2o_ice(ig,islope) > 0.) then
1183            albedo(ig,:,islope) = albedo_h2o_cap
1184            emis(ig,islope) = 1.
1185        endif
1186        ! CO2 ice
1187        if (co2_ice(ig,islope) > 0.) then
1188            albedo(ig,:,islope) = albedo_perennialco2(icap)
1189            emis(ig,islope) = emisice(icap)
1190        endif
1191        ! H2O frost
1192        if (qsurf(ig,igcm_h2o_ice,islope) > 0.) then
1193            albedo(ig,:,islope) = albedo_h2o_frost
1194            emis(ig,islope) = 1.
1195        endif
1196        ! CO2 frost
1197        if (qsurf(ig,igcm_co2,islope) > 0.) then
1198            albedo(ig,:,islope) = albedice(icap)
1199            emis(ig,islope) = emisice(icap)
1200        endif
1201    enddo
1202enddo
1203
[3770]1204! III_a.8 Orbital parameters update for start file
[3603]1205write(*,*) '> Setting the new orbital parameters'
[3498]1206if (evol_orbit_pem) call recomp_orb_param(i_myear,i_myear_leg)
[2779]1207
1208!------------------------
[3028]1209! III Output
[3317]1210!    III_b Write "restart.nc" and "restartfi.nc"
[3028]1211!------------------------
[3317]1212! III_b.1 Write "restart.nc"
[3042]1213ptimestep = iphysiq*daysec/real(day_step)/nsplit_phys ! dtphys/nsplit_phys
[3028]1214pday = day_ini
[3042]1215ztime_fin = time_phys
[2980]1216#ifndef CPP_1D
[3603]1217    write(*,*) '> Writing "restart.nc"'
[3442]1218    ! Correction on teta due to surface pressure changes
[3603]1219    allocate(pdyn(ip1jmp1))
1220    call gr_fi_dyn(1,ngrid,iip1,jjp1,ps_start0/ps_start,pdyn)
1221    do i = 1,ip1jmp1
[3779]1222        teta(i,:) = teta(i,:)*pdyn(i)**rcp
[3442]1223    enddo
1224    ! Correction on atmospheric pressure
[3603]1225    allocate(p(ip1jmp1,nlayer + 1))
1226    call gr_fi_dyn(1,ngrid,iip1,jjp1,ps_start,pdyn)
1227    call pression(ip1jmp1,ap,bp,pdyn,p)
[3442]1228    ! Correction on the mass of atmosphere
[3028]1229    call massdair(p,masse)
[3317]1230    call dynredem0("restart.nc",day_ini,phis)
[3628]1231    call dynredem1("restart.nc",time_0,vcov,ucov,teta,q,masse,pdyn)
[3603]1232    deallocate(ap,bp,p,pdyn)
[2980]1233#else
[3603]1234    write(*,*) '> Writing "restart1D.txt"'
[3571]1235    call writerestart1D('restart1D.txt',ps_start(1),tsurf(1,:),nlayer,size(tsurf,2),teta,ucov,vcov,nq,noms,qsurf(1,:,:),q)
[2980]1236#endif
[3591]1237deallocate(ps_start0,ps_start)
[2980]1238
[3317]1239! III_b.2 Write the "restartfi.nc"
[3603]1240write(*,*) '> Writing "restartfi.nc"'
[2842]1241#ifndef CPP_STD
[3317]1242    call physdem0("restartfi.nc",longitude,latitude,nsoilmx,ngrid, &
[3327]1243                  nlayer,nq,ptimestep,pday,0.,cell_area,albedodat, &
[3028]1244                  inertiedat,def_slope,subslope_dist)
[3327]1245    call physdem1("restartfi.nc",nsoilmx,ngrid,nlayer,nq,nqsoil,      &
[3114]1246                  ptimestep,ztime_fin,tsurf,tsoil,inertiesoil,        &
1247                  albedo,emis,q2,qsurf,qsoil,tauscaling,totcloudfrac, &
[3130]1248                  wstar,watercap,perennial_co2ice)
[2842]1249#else
[3591]1250    if (allocated(noms)) deallocate(noms)
1251    deallocate(qsurf,tsurf,tsoil,emis,watercap,watercaptag,albedo,inertiesoil)
[3317]1252    call physdem0("restartfi.nc",longitude,latitude,nsoilmx,ngrid, &
[3327]1253                  nlayer,nq,ptimestep,pday,time_phys,cell_area,    &
[3028]1254                  albedo_bareground,inertiedat,zmea,zstd,zsig,zgam,zthe)
[3327]1255    call physdem1("restartfi.nc",nsoilmx,ngrid,nlayer,nq,nqsoil,       &
[3114]1256                  ptimestep,ztime_fin,tsurf,tsoil,emis,q2,qsurf,qsoil, &
1257                  cloudfrac,totcloudfrac,hice,rnat,pctsrf_sic,tslab,   &
1258                  tsea_ice,sea_ice)
[2842]1259#endif
1260
[2794]1261!------------------------
1262! III Output
[3161]1263!    III_c Write the "restartpem.nc"
[2794]1264!------------------------
[3603]1265write(*,*) '> Writing "restartpem.nc"'
[3778]1266if (layering_algo) nb_str_max = get_nb_str_max(layerings_map,ngrid,nslope) ! Get the maximum number of "stratum" in the layerings_mapication (layerings)
[3206]1267call pemdem0("restartpem.nc",longitude,latitude,cell_area,ngrid,nslope,def_slope,subslope_dist)
[3537]1268call pemdem1("restartpem.nc",i_myear,nsoilmx_PEM,ngrid,nslope,tsoil_PEM,TI_PEM,icetable_depth,icetable_thickness,ice_porefilling, &
[3778]1269             co2_adsorbed_phys,h2o_adsorbed_phys,h2o_ice,layerings_map)
[2779]1270
[3498]1271call info_PEM(i_myear_leg,stopPEM,i_myear,n_myear)
[3149]1272
[3616]1273write(*,*)
[3738]1274write(*,*) '****** PEM final information *******'
[3770]1275write(*,'(a,f16.4,a)') " + The PEM leg has run for ", i_myear_leg, " Martian years."
1276write(*,'(a,f16.4,a,f16.4,a)') " + The chained simulation has run for ", i_myear, " Martian years =", i_myear*convert_years, " Earth years."
1277write(*,'(a,f16.4,a)') " + The reached date is now ", (year_bp_ini + i_myear)*convert_years, " Earth years."
[3616]1278write(*,*) "+ PEM: so far, so good!"
[3673]1279write(*,*) '************************************'
[2794]1280
[3319]1281if (layering_algo) then
1282    do islope = 1,nslope
1283        do i = 1,ngrid
[3778]1284            call del_layering(layerings_map(i,islope))
[3319]1285        enddo
[3297]1286    enddo
[3319]1287endif
[3591]1288deallocate(q,longitude,latitude,cell_area,tsoil_PEM)
[3778]1289deallocate(co2_ice,h2o_ice,layerings_map)
[3028]1290!----------------------------- END OUTPUT ------------------------------
[2897]1291
[2779]1292END PROGRAM pem
Note: See TracBrowser for help on using the repository browser.