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