[3149] | 1 | MODULE glaciers_mod |
---|
[2995] | 2 | |
---|
[3149] | 3 | implicit none |
---|
| 4 | |
---|
[3161] | 5 | ! Flags for ice management |
---|
[3308] | 6 | logical :: h2oice_flow ! True by default, to compute H2O ice flow. Read in "run_PEM.def" |
---|
| 7 | logical :: co2ice_flow ! True by default, to compute CO2 ice flow. Read in "run_PEM.def" |
---|
[3161] | 8 | logical :: metam_h2oice ! False by default, to compute H2O ice metamorphism. Read in "run_PEM.def" |
---|
| 9 | logical :: metam_co2ice ! False by default, to compute CO2 ice metamorphism. Read in "run_PEM.def" |
---|
| 10 | |
---|
| 11 | ! Thresholds for ice management |
---|
[3256] | 12 | real, save :: inf_h2oice_threshold ! To consider the amount of H2O ice as an infinite reservoir [kg.m-2] |
---|
| 13 | real, save :: metam_h2oice_threshold ! To consider frost is becoming perennial H2O ice [kg.m-2] |
---|
| 14 | real, save :: metam_co2ice_threshold ! To consider frost is becoming perennial CO2 ice [kg.m-2] |
---|
[3161] | 15 | |
---|
[3149] | 16 | !======================================================================= |
---|
[2995] | 17 | contains |
---|
[3149] | 18 | !======================================================================= |
---|
[2995] | 19 | |
---|
[3149] | 20 | SUBROUTINE flow_co2glaciers(timelen,ngrid,nslope,iflat,subslope_dist,def_slope_mean,vmr_co2_PEM,ps_PCM,global_avg_ps_PCM,global_avg_ps_PEM,co2ice,flag_co2flow,flag_co2flow_mesh) |
---|
[2995] | 21 | |
---|
| 22 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 23 | !!! |
---|
| 24 | !!! Purpose: Main for CO2 glaciers evolution: compute maximum thickness, and do |
---|
| 25 | !!! the ice transfer |
---|
[3532] | 26 | !!! |
---|
| 27 | !!! |
---|
[2995] | 28 | !!! Author: LL |
---|
| 29 | !!! |
---|
| 30 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 31 | |
---|
[3149] | 32 | implicit none |
---|
[2995] | 33 | |
---|
[3527] | 34 | ! Inputs: |
---|
| 35 | !-------- |
---|
| 36 | integer, intent(in) :: timelen, ngrid, nslope, iflat ! number of time sample, physical points, subslopes, index of the flat subslope |
---|
| 37 | real, dimension(ngrid,nslope), intent(in) :: subslope_dist ! Physical points x Slopes: Distribution of the subgrid slopes |
---|
| 38 | real, dimension(ngrid), intent(in) :: def_slope_mean ! Physical points: values of the sub grid slope angles |
---|
| 39 | real, dimension(ngrid,timelen), intent(in) :: vmr_co2_PEM ! Physical x Time field : VMR of co2 in the first layer [mol/mol] |
---|
| 40 | real, dimension(ngrid,timelen), intent(in) :: ps_PCM ! Physical x Time field: surface pressure given by the PCM [Pa] |
---|
| 41 | real, intent(in) :: global_avg_ps_PCM ! Global averaged surface pressure from the PCM [Pa] |
---|
| 42 | real, intent(in) :: global_avg_ps_PEM ! global averaged surface pressure during the PEM iteration [Pa] |
---|
[2995] | 43 | |
---|
| 44 | ! Ouputs: |
---|
[3527] | 45 | !-------- |
---|
| 46 | real, dimension(ngrid,nslope), intent(inout) :: co2ice ! Physical x Slope field: co2 ice on the subgrid slopes [kg/m^2] |
---|
| 47 | real, dimension(ngrid,nslope), intent(inout) :: flag_co2flow ! flag to see if there is flow on the subgrid slopes |
---|
| 48 | real, dimension(ngrid), intent(inout) :: flag_co2flow_mesh ! same but within the mesh |
---|
[3532] | 49 | ! Local |
---|
[3527] | 50 | !------ |
---|
| 51 | real, dimension(ngrid,nslope) :: Tcond ! Physical field: CO2 condensation temperature [K] |
---|
[3532] | 52 | real, dimension(ngrid,nslope) :: hmax ! Physical x Slope field: maximum thickness for co2 glacier before flow |
---|
[2995] | 53 | |
---|
[3345] | 54 | write(*,*) "Flow of CO2 glaciers" |
---|
[3149] | 55 | call computeTcondCO2(timelen,ngrid,nslope,vmr_co2_PEM,ps_PCM,global_avg_ps_PCM,global_avg_ps_PEM,Tcond) |
---|
| 56 | call compute_hmaxglaciers(ngrid,nslope,iflat,def_slope_mean,Tcond,"co2",hmax) |
---|
| 57 | call transfer_ice_duringflow(ngrid,nslope,iflat, subslope_dist,def_slope_mean,hmax,Tcond,"co2",co2ice,flag_co2flow,flag_co2flow_mesh) |
---|
[2995] | 58 | |
---|
[3149] | 59 | END SUBROUTINE flow_co2glaciers |
---|
[2995] | 60 | |
---|
[3149] | 61 | !======================================================================= |
---|
[2995] | 62 | |
---|
[3149] | 63 | SUBROUTINE flow_h2oglaciers(timelen,ngrid,nslope,iflat,subslope_dist,def_slope_mean,Tice,h2oice,flag_h2oflow,flag_h2oflow_mesh) |
---|
[2995] | 64 | |
---|
| 65 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 66 | !!! |
---|
| 67 | !!! Purpose: Main for H2O glaciers evolution: compute maximum thickness, and do |
---|
| 68 | !!! the ice transfer |
---|
[3532] | 69 | !!! |
---|
| 70 | !!! |
---|
[2995] | 71 | !!! Author: LL |
---|
| 72 | !!! |
---|
| 73 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 74 | |
---|
[3149] | 75 | implicit none |
---|
[2995] | 76 | |
---|
| 77 | ! arguments |
---|
| 78 | ! --------- |
---|
| 79 | |
---|
| 80 | ! Inputs: |
---|
[3527] | 81 | integer, intent(in) :: timelen, ngrid, nslope, iflat ! number of time sample, physical points, subslopes, index of the flat subslope |
---|
| 82 | real, dimension(ngrid,nslope), intent(in) :: subslope_dist ! Physical points x Slopes : Distribution of the subgrid slopes |
---|
| 83 | real, dimension(ngrid), intent(in) :: def_slope_mean ! Slopes: values of the sub grid slope angles |
---|
| 84 | real, dimension(ngrid,nslope), intent(in) :: Tice ! Ice Temperature [K] |
---|
[2995] | 85 | ! Ouputs: |
---|
[3527] | 86 | real, dimension(ngrid,nslope), intent(inout) :: h2oice ! Physical x Slope field: co2 ice on the subgrid slopes [kg/m^2] |
---|
| 87 | real, dimension(ngrid,nslope), intent(inout) :: flag_h2oflow ! flag to see if there is flow on the subgrid slopes |
---|
| 88 | real, dimension(ngrid), intent(inout) :: flag_h2oflow_mesh ! same but within the mesh |
---|
[3532] | 89 | ! Local |
---|
| 90 | real, dimension(ngrid,nslope) :: hmax ! Physical x Slope field: maximum thickness for co2 glacier before flow |
---|
[2995] | 91 | |
---|
[3149] | 92 | write(*,*) "Flow of H2O glaciers" |
---|
| 93 | call compute_hmaxglaciers(ngrid,nslope,iflat,def_slope_mean,Tice,"h2o",hmax) |
---|
| 94 | call transfer_ice_duringflow(ngrid,nslope,iflat, subslope_dist,def_slope_mean,hmax,Tice,"h2o",h2oice,flag_h2oflow,flag_h2oflow_mesh) |
---|
[2995] | 95 | |
---|
[3149] | 96 | END SUBROUTINE flow_h2oglaciers |
---|
[2995] | 97 | |
---|
[3149] | 98 | !======================================================================= |
---|
[2995] | 99 | |
---|
[3149] | 100 | SUBROUTINE compute_hmaxglaciers(ngrid,nslope,iflat,def_slope_mean,Tice,name_ice,hmax) |
---|
[2995] | 101 | |
---|
[3527] | 102 | use ice_table_mod, only: rho_ice |
---|
[3076] | 103 | use abort_pem_mod, only: abort_pem |
---|
[3082] | 104 | #ifndef CPP_STD |
---|
| 105 | use comcstfi_h, only: pi, g |
---|
| 106 | #else |
---|
| 107 | use comcstfi_mod, only: pi, g |
---|
| 108 | #endif |
---|
[2995] | 109 | |
---|
| 110 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 111 | !!! |
---|
[3527] | 112 | !!! Purpose: Compute the maximum thickness of CO2 and H2O glaciers given a slope angle before initating flow |
---|
[3532] | 113 | !!! |
---|
[2995] | 114 | !!! Author: LL,based on work by A.Grau Galofre (LPG) and Isaac Smith (JGR Planets 2022) |
---|
| 115 | !!! |
---|
| 116 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 117 | |
---|
[3149] | 118 | implicit none |
---|
[2995] | 119 | |
---|
| 120 | ! Inputs |
---|
[3527] | 121 | ! ------ |
---|
| 122 | integer, intent(in) :: ngrid, nslope ! # of grid points and subslopes |
---|
| 123 | integer, intent(in) :: iflat ! index of the flat subslope |
---|
| 124 | real, dimension(nslope), intent(in) :: def_slope_mean ! Slope field: Values of the subgrid slope angles [deg] |
---|
| 125 | real, dimension(ngrid,nslope), intent(in) :: Tice ! Physical field: ice temperature [K] |
---|
| 126 | character(3), intent(in) :: name_ice ! Nature of ice |
---|
[2995] | 127 | ! Outputs |
---|
[3527] | 128 | ! ------- |
---|
| 129 | real, dimension(ngrid,nslope), intent(out) :: hmax ! Physical grid x Slope field: maximum thickness before flaw [m] |
---|
[2995] | 130 | ! Local |
---|
[3527] | 131 | ! ----- |
---|
| 132 | real :: tau_d ! characteristic basal drag, understood as the stress that an ice mass flowing under its weight balanced by viscosity. Value obtained from I.Smith |
---|
| 133 | integer :: ig, islope ! loop variables |
---|
| 134 | real :: slo_angle |
---|
[2995] | 135 | |
---|
[3527] | 136 | select case (trim(adjustl(name_ice))) |
---|
| 137 | case('h2o') |
---|
| 138 | tau_d = 1.e5 |
---|
| 139 | case('co2') |
---|
[2995] | 140 | tau_d = 5.e3 |
---|
[3527] | 141 | case default |
---|
| 142 | call abort_pem("compute_hmaxglaciers","Type of ice not known!",1) |
---|
| 143 | end select |
---|
[2995] | 144 | |
---|
[3527] | 145 | do ig = 1,ngrid |
---|
| 146 | do islope = 1,nslope |
---|
| 147 | if (islope == iflat) then |
---|
[2995] | 148 | hmax(ig,islope) = 1.e8 |
---|
[3527] | 149 | else |
---|
[2995] | 150 | slo_angle = abs(def_slope_mean(islope)*pi/180.) |
---|
[3527] | 151 | hmax(ig,islope) = tau_d/(rho_ice(Tice(ig,islope),name_ice)*g*slo_angle) |
---|
| 152 | endif |
---|
| 153 | enddo |
---|
| 154 | enddo |
---|
| 155 | |
---|
[3149] | 156 | END SUBROUTINE compute_hmaxglaciers |
---|
[2995] | 157 | |
---|
[3149] | 158 | !======================================================================= |
---|
[2995] | 159 | |
---|
[3149] | 160 | SUBROUTINE transfer_ice_duringflow(ngrid,nslope,iflat,subslope_dist,def_slope_mean,hmax,Tice,name_ice,qice,flag_flow,flag_flowmesh) |
---|
[2995] | 161 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 162 | !!! |
---|
| 163 | !!! Purpose: Transfer the excess of ice from one subslope to another |
---|
| 164 | !!! No transfer between mesh at the time |
---|
| 165 | !!! Author: LL |
---|
| 166 | !!! |
---|
| 167 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 168 | |
---|
[3527] | 169 | use ice_table_mod, only: rho_ice |
---|
[3076] | 170 | use abort_pem_mod, only: abort_pem |
---|
[3082] | 171 | #ifndef CPP_STD |
---|
| 172 | use comcstfi_h, only: pi |
---|
| 173 | #else |
---|
| 174 | use comcstfi_mod, only: pi |
---|
| 175 | #endif |
---|
[2995] | 176 | |
---|
| 177 | implicit none |
---|
| 178 | |
---|
| 179 | ! Inputs |
---|
[3527] | 180 | !------- |
---|
| 181 | integer, intent(in) :: ngrid, nslope ! # of physical points and subslope |
---|
| 182 | integer, intent(in) :: iflat ! index of the flat subslope |
---|
| 183 | real, dimension(ngrid,nslope), intent(in) :: subslope_dist ! Distribution of the subgrid slopes within the mesh |
---|
| 184 | real, dimension(nslope), intent(in) :: def_slope_mean ! values of the subgrid slopes |
---|
| 185 | real, dimension(ngrid,nslope), intent(in) :: hmax ! maximum height of the glaciers before initiating flow [m] |
---|
| 186 | real, dimension(ngrid,nslope), intent(in) :: Tice ! Ice temperature[K] |
---|
| 187 | character(3), intent(in) :: name_ice ! Nature of the ice |
---|
[2995] | 188 | ! Outputs |
---|
[3527] | 189 | !-------- |
---|
| 190 | real, dimension(ngrid,nslope), intent(inout) :: qice ! CO2 in the subslope [kg/m^2] |
---|
| 191 | real, dimension(ngrid,nslope), intent(inout) :: flag_flow ! boolean to check if there is flow on a subgrid slope |
---|
| 192 | real, dimension(ngrid), intent(inout) :: flag_flowmesh ! boolean to check if there is flow in the mesh |
---|
[2995] | 193 | ! Local |
---|
[3527] | 194 | !------ |
---|
| 195 | integer :: ig, islope ! loop |
---|
| 196 | integer :: iaval ! ice will be transfered here |
---|
[2995] | 197 | |
---|
[3527] | 198 | do ig = 1,ngrid |
---|
| 199 | do islope = 1,nslope |
---|
| 200 | if (islope /= iflat) then ! ice can be infinite on flat ground |
---|
[2995] | 201 | ! First: check that CO2 ice must flow (excess of ice on the slope), ice can accumulate infinitely on flat ground |
---|
[3527] | 202 | if (qice(ig,islope) >= rho_ice(Tice(ig,islope),'h2o')*hmax(ig,islope)*cos(pi*def_slope_mean(islope)/180.)) then |
---|
[2995] | 203 | ! Second: determine the flatest slopes possible: |
---|
[3527] | 204 | if (islope > iflat) then |
---|
| 205 | iaval=islope-1 |
---|
| 206 | else |
---|
| 207 | iaval = islope + 1 |
---|
| 208 | endif |
---|
| 209 | do while (iaval /= iflat .and. subslope_dist(ig,iaval) == 0) |
---|
| 210 | if (iaval > iflat) then |
---|
| 211 | iaval = iaval - 1 |
---|
| 212 | else |
---|
| 213 | iaval = iaval + 1 |
---|
| 214 | endif |
---|
[2995] | 215 | enddo |
---|
[3527] | 216 | qice(ig,iaval) = qice(ig,iaval) + (qice(ig,islope) - rho_ice(Tice(ig,islope),'h2o')*hmax(ig,islope)*cos(pi*def_slope_mean(islope)/180.)) & |
---|
| 217 | *subslope_dist(ig,islope)/subslope_dist(ig,iaval)*cos(pi*def_slope_mean(iaval)/180.)/cos(pi*def_slope_mean(islope)/180.) |
---|
[2995] | 218 | |
---|
[3527] | 219 | qice(ig,islope) = rho_ice(Tice(ig,islope),'h2o')*hmax(ig,islope)*cos(pi*def_slope_mean(islope)/180.) |
---|
| 220 | flag_flow(ig,islope) = 1. |
---|
| 221 | flag_flowmesh(ig) = 1. |
---|
| 222 | endif ! co2ice > hmax |
---|
| 223 | endif ! iflat |
---|
[3532] | 224 | enddo !islope |
---|
[3527] | 225 | enddo !ig |
---|
[2995] | 226 | |
---|
[3149] | 227 | END SUBROUTINE |
---|
[2995] | 228 | |
---|
[3149] | 229 | !======================================================================= |
---|
| 230 | |
---|
| 231 | SUBROUTINE computeTcondCO2(timelen,ngrid,nslope,vmr_co2_PEM,ps_PCM,global_avg_ps_PCM,global_avg_ps_PEM,Tcond) |
---|
[2995] | 232 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
| 233 | !!! |
---|
| 234 | !!! Purpose: Compute CO2 condensation temperature |
---|
| 235 | !!! |
---|
| 236 | !!! Author: LL |
---|
| 237 | !!! |
---|
| 238 | !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! |
---|
[3532] | 239 | |
---|
[3527] | 240 | use constants_marspem_mod, only: alpha_clap_co2,beta_clap_co2 |
---|
[2995] | 241 | |
---|
[3532] | 242 | implicit none |
---|
[2995] | 243 | |
---|
| 244 | ! arguments: |
---|
| 245 | ! ---------- |
---|
| 246 | |
---|
| 247 | ! INPUT |
---|
[3149] | 248 | integer, intent(in) :: timelen, ngrid, nslope ! # of timesample, physical points, subslopes |
---|
| 249 | real, dimension(ngrid,timelen), intent(in) :: vmr_co2_PEM ! Physical points x times field: VMR of CO2 in the first layer [mol/mol] |
---|
| 250 | real, dimension(ngrid,timelen), intent(in) :: ps_PCM ! Physical points x times field: surface pressure in the PCM [Pa] |
---|
| 251 | real, intent(in) :: global_avg_ps_PCM ! Global averaged surfacepressure in the PCM [Pa] |
---|
| 252 | real, intent(in) :: global_avg_ps_PEM ! Global averaged surface pressure computed during the PEM iteration |
---|
| 253 | |
---|
[2995] | 254 | ! OUTPUT |
---|
[3532] | 255 | real, dimension(ngrid,nslope), intent(out) :: Tcond ! Physical points: condensation temperature of CO2, yearly averaged |
---|
[2995] | 256 | |
---|
| 257 | ! LOCAL |
---|
[3149] | 258 | integer :: ig, it ! For loop |
---|
| 259 | real :: ave ! Intermediate to compute average |
---|
[2995] | 260 | |
---|
[3149] | 261 | do ig = 1,ngrid |
---|
| 262 | ave = 0 |
---|
| 263 | do it = 1,timelen |
---|
| 264 | ave = ave + beta_clap_co2/(alpha_clap_co2 - log(vmr_co2_PEM(ig,it)*ps_PCM(ig,it)*global_avg_ps_PCM/global_avg_ps_PEM/100)) |
---|
| 265 | enddo |
---|
| 266 | Tcond(ig,:) = ave/timelen |
---|
| 267 | enddo |
---|
[2995] | 268 | |
---|
[3149] | 269 | END SUBROUTINE computeTcondCO2 |
---|
[2995] | 270 | |
---|
[3149] | 271 | END MODULE glaciers_mod |
---|