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