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

Last change on this file since 3578 was 3578, checked in by jbclement, 6 days ago

COMMON:
Follow-up of r3576:

  • Status result is added to the file "version_stat-diff.txt", previously called "version_diff.txt";
  • Addition of the functionnality for the program 'Reshape' dedicated to the PEM;
  • Improvements of the visualization format.

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