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

Last change on this file since 3386 was 3386, checked in by jbclement, 21 months ago

PEM:

  • Correction for the initialization in 1D due to the change of 'ecritphy' into 'outputs_per_sol' in r3369.
  • Correction in the condition to check the normal end of a PCM/PEM run in "jobP*M.slurm" files.
  • Correction of files cleaning in the launching script.
  • Improvement of the launching script to take into account the modification of parameters for a relaunch.

JBC

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