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

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

PEM:

  • New way to manage the pressure: now the PEM manages only the average pressure and keeps the pressure deviation with the instantaneous pressure from the start to reconstruct the pressure at the end ('ps_avg = ps_start + ps_dev'). As a consequence, everything related to pressure in the PEM is modified accordingly.
  • Surface temperatures management is now simpler. It follows the strategy for the pressure (and soil temperature) described above.
  • Soil temperatures are now adapted to match the surface temperature changes occured during the PEM by modifying the soil temperature deviation at the end.
  • Few simplifications/optimizations: notably, the two PCM years are now read in one go in 'read_data_PCM_mod.F90' and only the needed variables are extracted.
  • Deletion of unused variables and unnecessary intermediate variables (memory saving and loop deletion in some cases).
  • Renaming of variables and subroutines to make everything clearer. In particular, the suffixes: '_avg' = average, '_start' = PCM start file, '_dev' = deviation, '_ini' or '0' = initial, '_dyn' = dynamical grid, '_timeseries' = daily average of last PCM year.
  • Cosmetic cleanings for readability.

JBC

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