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

Last change on this file since 3599 was 3599, checked in by jbclement, 29 hours ago

PEM:
Few small optimizations.
JBC

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