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