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

Last change on this file since 3785 was 3785, checked in by jbclement, 3 weeks ago

PEM:

  • Rework of ice sublimation in the layering algorithm to sublimate ice and to lift dust at the same time
  • Addition of the possibility to make H2O ice or mixed dust/H2O ice lag when CO2 ice is sublimating
  • Addition of a lag resistance for CO2 ice sublimation
  • Consideration of a 'wet_lag_porosity' for H2O ice lag layer
  • Consideration of a threshold under which H2O ice lag layer is considered patchy
  • Deletion of unused function 'thickness_toplag'

JBC

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