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

Last change on this file since 3490 was 3490, checked in by jbclement, 2 weeks ago

PEM:
Correction of initialization of 'icetable_dynamic' and the management of the flag.
JBC

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