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

Last change on this file since 3410 was 3403, checked in by jbclement, 11 months ago

PEM:
Addition in the launching script of the possibility to submit a job with PBS/TORQUE + Modification to make the time limit detection in "pem.F90" work with PBS/TORQUE + Update of the headers of .job files.
JBC

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