| 1 | MODULE surf_ice |
|---|
| 2 | !----------------------------------------------------------------------- |
|---|
| 3 | ! NAME |
|---|
| 4 | ! surf_ice |
|---|
| 5 | ! |
|---|
| 6 | ! DESCRIPTION |
|---|
| 7 | ! Surface ice management. |
|---|
| 8 | ! |
|---|
| 9 | ! AUTHORS & DATE |
|---|
| 10 | ! JB Clement, 12/2025 |
|---|
| 11 | ! |
|---|
| 12 | ! NOTES |
|---|
| 13 | ! |
|---|
| 14 | !----------------------------------------------------------------------- |
|---|
| 15 | |
|---|
| 16 | ! DEPENDENCIES |
|---|
| 17 | ! ------------ |
|---|
| 18 | use numerics, only: dp, qp, di, k4 |
|---|
| 19 | |
|---|
| 20 | ! DECLARATION |
|---|
| 21 | ! ----------- |
|---|
| 22 | implicit none |
|---|
| 23 | |
|---|
| 24 | ! PARAMETERS |
|---|
| 25 | ! ---------- |
|---|
| 26 | real(dp), parameter :: rho_co2ice = 1650._dp ! Density of CO2 ice [kg.m-3] |
|---|
| 27 | real(dp), parameter :: rho_h2oice = 920._dp ! Density of H2O ice [kg.m-3] |
|---|
| 28 | real(dp), protected :: threshold_h2oice_cap ! Threshold to consider H2O ice as infinite reservoir [kg.m-2] |
|---|
| 29 | real(dp), protected :: h2oice_huge_ini ! Initial value for huge reservoirs of H2O ice [kg/m^2] |
|---|
| 30 | logical(k4), dimension(:), allocatable, protected :: is_h2o_perice_PCM ! Flag for the presence of perennial H2O ice in the PCM at the beginning [kg.m-2] |
|---|
| 31 | real(dp), dimension(:,:), allocatable, protected :: co2_perice_PCM ! Perennial CO2 ice in the PCM at the beginning [kg.m-2] |
|---|
| 32 | |
|---|
| 33 | contains |
|---|
| 34 | !+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ |
|---|
| 35 | |
|---|
| 36 | !======================================================================= |
|---|
| 37 | SUBROUTINE ini_surf_ice() |
|---|
| 38 | !----------------------------------------------------------------------- |
|---|
| 39 | ! NAME |
|---|
| 40 | ! ini_surf_ice |
|---|
| 41 | ! |
|---|
| 42 | ! DESCRIPTION |
|---|
| 43 | ! Initialize surface ice arrays. |
|---|
| 44 | ! |
|---|
| 45 | ! AUTHORS & DATE |
|---|
| 46 | ! JB Clement 12/2025 |
|---|
| 47 | ! |
|---|
| 48 | ! NOTES |
|---|
| 49 | ! |
|---|
| 50 | !----------------------------------------------------------------------- |
|---|
| 51 | |
|---|
| 52 | ! DEPENDENCIES |
|---|
| 53 | ! ------------ |
|---|
| 54 | use geometry, only: ngrid, nslope |
|---|
| 55 | |
|---|
| 56 | ! DECLARATION |
|---|
| 57 | ! ----------- |
|---|
| 58 | implicit none |
|---|
| 59 | |
|---|
| 60 | ! CODE |
|---|
| 61 | ! ---- |
|---|
| 62 | if (.not. allocated(is_h2o_perice_PCM)) allocate(is_h2o_perice_PCM(ngrid)) |
|---|
| 63 | if (.not. allocated(co2_perice_PCM)) allocate(co2_perice_PCM(ngrid,nslope)) |
|---|
| 64 | is_h2o_perice_PCM(:) = .false. |
|---|
| 65 | co2_perice_PCM(:,:) = 0._dp |
|---|
| 66 | |
|---|
| 67 | END SUBROUTINE ini_surf_ice |
|---|
| 68 | !======================================================================= |
|---|
| 69 | |
|---|
| 70 | !======================================================================= |
|---|
| 71 | SUBROUTINE end_surf_ice() |
|---|
| 72 | !----------------------------------------------------------------------- |
|---|
| 73 | ! NAME |
|---|
| 74 | ! end_surf_ice |
|---|
| 75 | ! |
|---|
| 76 | ! DESCRIPTION |
|---|
| 77 | ! Deallocate surface ice arrays. |
|---|
| 78 | ! |
|---|
| 79 | ! AUTHORS & DATE |
|---|
| 80 | ! JB Clement, 12/2025 |
|---|
| 81 | ! |
|---|
| 82 | ! NOTES |
|---|
| 83 | ! |
|---|
| 84 | !----------------------------------------------------------------------- |
|---|
| 85 | |
|---|
| 86 | ! DECLARATION |
|---|
| 87 | ! ----------- |
|---|
| 88 | implicit none |
|---|
| 89 | |
|---|
| 90 | ! CODE |
|---|
| 91 | ! ---- |
|---|
| 92 | if (allocated(is_h2o_perice_PCM)) deallocate(is_h2o_perice_PCM) |
|---|
| 93 | if (allocated(co2_perice_PCM)) deallocate(co2_perice_PCM) |
|---|
| 94 | |
|---|
| 95 | END SUBROUTINE end_surf_ice |
|---|
| 96 | !======================================================================= |
|---|
| 97 | |
|---|
| 98 | !======================================================================= |
|---|
| 99 | SUBROUTINE set_surf_ice_config(threshold_h2oice_cap_in,h2oice_huge_ini_in) |
|---|
| 100 | !----------------------------------------------------------------------- |
|---|
| 101 | ! NAME |
|---|
| 102 | ! set_surf_ice_config |
|---|
| 103 | ! |
|---|
| 104 | ! DESCRIPTION |
|---|
| 105 | ! Setter for 'surf_ice' configuration parameters. |
|---|
| 106 | ! |
|---|
| 107 | ! AUTHORS & DATE |
|---|
| 108 | ! JB Clement, 02/2026 |
|---|
| 109 | ! |
|---|
| 110 | ! NOTES |
|---|
| 111 | ! |
|---|
| 112 | !----------------------------------------------------------------------- |
|---|
| 113 | |
|---|
| 114 | ! DEPENDENCIES |
|---|
| 115 | ! ------------ |
|---|
| 116 | use utility, only: real2str |
|---|
| 117 | use display, only: print_msg |
|---|
| 118 | |
|---|
| 119 | ! DECLARATION |
|---|
| 120 | ! ----------- |
|---|
| 121 | implicit none |
|---|
| 122 | |
|---|
| 123 | ! ARGUMENTS |
|---|
| 124 | ! --------- |
|---|
| 125 | real(dp), intent(in) :: threshold_h2oice_cap_in, h2oice_huge_ini_in |
|---|
| 126 | |
|---|
| 127 | ! CODE |
|---|
| 128 | ! ---- |
|---|
| 129 | threshold_h2oice_cap = threshold_h2oice_cap_in |
|---|
| 130 | h2oice_huge_ini = h2oice_huge_ini_in |
|---|
| 131 | call print_msg('threshold_h2oice_cap = '//real2str(threshold_h2oice_cap)) |
|---|
| 132 | call print_msg('h2oice_huge_ini = '//real2str(h2oice_huge_ini)) |
|---|
| 133 | |
|---|
| 134 | END SUBROUTINE set_surf_ice_config |
|---|
| 135 | !======================================================================= |
|---|
| 136 | |
|---|
| 137 | !======================================================================= |
|---|
| 138 | SUBROUTINE set_co2_perice_PCM(co2_perice_PCM_in) |
|---|
| 139 | !----------------------------------------------------------------------- |
|---|
| 140 | ! NAME |
|---|
| 141 | ! set_co2_perice_PCM |
|---|
| 142 | ! |
|---|
| 143 | ! DESCRIPTION |
|---|
| 144 | ! Setter for 'co2_perice_PCM'. |
|---|
| 145 | ! |
|---|
| 146 | ! AUTHORS & DATE |
|---|
| 147 | ! JB Clement, 12/2025 |
|---|
| 148 | ! |
|---|
| 149 | ! NOTES |
|---|
| 150 | ! |
|---|
| 151 | !----------------------------------------------------------------------- |
|---|
| 152 | |
|---|
| 153 | ! DECLARATION |
|---|
| 154 | ! ----------- |
|---|
| 155 | implicit none |
|---|
| 156 | |
|---|
| 157 | ! ARGUMENTS |
|---|
| 158 | ! --------- |
|---|
| 159 | real(dp), dimension(:,:), intent(in) :: co2_perice_PCM_in |
|---|
| 160 | |
|---|
| 161 | ! CODE |
|---|
| 162 | ! ---- |
|---|
| 163 | co2_perice_PCM(:,:) = co2_perice_PCM_in(:,:) |
|---|
| 164 | |
|---|
| 165 | END SUBROUTINE set_co2_perice_PCM |
|---|
| 166 | !======================================================================= |
|---|
| 167 | |
|---|
| 168 | !======================================================================= |
|---|
| 169 | SUBROUTINE set_is_h2o_perice_PCM(is_h2o_perice_PCM_in) |
|---|
| 170 | !----------------------------------------------------------------------- |
|---|
| 171 | ! NAME |
|---|
| 172 | ! set_is_h2o_perice_PCM |
|---|
| 173 | ! |
|---|
| 174 | ! DESCRIPTION |
|---|
| 175 | ! Setter for 'is_h2o_perice_PCM'. |
|---|
| 176 | ! |
|---|
| 177 | ! AUTHORS & DATE |
|---|
| 178 | ! JB Clement, 12/2025 |
|---|
| 179 | ! |
|---|
| 180 | ! NOTES |
|---|
| 181 | ! |
|---|
| 182 | !----------------------------------------------------------------------- |
|---|
| 183 | |
|---|
| 184 | ! DECLARATION |
|---|
| 185 | ! ----------- |
|---|
| 186 | implicit none |
|---|
| 187 | |
|---|
| 188 | ! ARGUMENTS |
|---|
| 189 | ! --------- |
|---|
| 190 | real(dp), dimension(:), intent(in) :: is_h2o_perice_PCM_in |
|---|
| 191 | |
|---|
| 192 | ! CODE |
|---|
| 193 | ! ---- |
|---|
| 194 | where(is_h2o_perice_PCM_in > 0.5_dp) |
|---|
| 195 | is_h2o_perice_PCM = .true. |
|---|
| 196 | else where |
|---|
| 197 | is_h2o_perice_PCM = .false. |
|---|
| 198 | end where |
|---|
| 199 | |
|---|
| 200 | END SUBROUTINE set_is_h2o_perice_PCM |
|---|
| 201 | !======================================================================= |
|---|
| 202 | |
|---|
| 203 | !======================================================================= |
|---|
| 204 | SUBROUTINE build4PCM_perice(h2o_ice,co2_ice,is_h2o_perice,h2o_ice4PCM,co2_ice4PCM) |
|---|
| 205 | !----------------------------------------------------------------------- |
|---|
| 206 | ! NAME |
|---|
| 207 | ! build4PCM_perice |
|---|
| 208 | ! |
|---|
| 209 | ! DESCRIPTION |
|---|
| 210 | ! Reconstructs perennial ice and frost for the PCM. |
|---|
| 211 | ! |
|---|
| 212 | ! AUTHORS & DATE |
|---|
| 213 | ! JB Clement, 12/2025 |
|---|
| 214 | ! |
|---|
| 215 | ! NOTES |
|---|
| 216 | ! |
|---|
| 217 | !----------------------------------------------------------------------- |
|---|
| 218 | |
|---|
| 219 | ! DEPENDENCIES |
|---|
| 220 | ! ------------ |
|---|
| 221 | use geometry, only: ngrid |
|---|
| 222 | use frost, only: h2o_frost4PCM |
|---|
| 223 | use slopes, only: subslope_dist, def_slope_mean |
|---|
| 224 | use maths, only: pi |
|---|
| 225 | use display, only: print_msg |
|---|
| 226 | |
|---|
| 227 | ! DECLARATION |
|---|
| 228 | ! ----------- |
|---|
| 229 | implicit none |
|---|
| 230 | |
|---|
| 231 | ! ARGUMENTS |
|---|
| 232 | ! --------- |
|---|
| 233 | real(dp), dimension(:,:), intent(inout) :: h2o_ice, co2_ice |
|---|
| 234 | logical(k4), dimension(:), intent(out) :: is_h2o_perice ! H2O perennial ice flag |
|---|
| 235 | real(dp), dimension(:,:), intent(out) :: h2o_ice4PCM, co2_ice4PCM ! Ice for PCM |
|---|
| 236 | |
|---|
| 237 | ! LOCAL VARIABLES |
|---|
| 238 | ! --------------- |
|---|
| 239 | integer(di) :: i |
|---|
| 240 | |
|---|
| 241 | ! CODE |
|---|
| 242 | ! ---- |
|---|
| 243 | call print_msg('> Building surface ice/frost for the PCM') |
|---|
| 244 | co2_ice4PCM(:,:) = co2_ice(:,:) |
|---|
| 245 | h2o_ice4PCM(:,:) = 0._dp ! Because in the Mars PCM, only the variation of perennial H2O ice is monitored, not the absolute quantity |
|---|
| 246 | do i = 1,ngrid |
|---|
| 247 | ! Is H2O ice still considered as an infinite reservoir for the PCM? |
|---|
| 248 | if (sum(h2o_ice(i,:)*subslope_dist(i,:)/cos(pi*def_slope_mean(:)/180._dp)) > threshold_h2oice_cap) then |
|---|
| 249 | ! There is enough ice to be considered as an infinite reservoir |
|---|
| 250 | is_h2o_perice(i) = .true. |
|---|
| 251 | else |
|---|
| 252 | ! Too little ice to be considered as an infinite reservoir so ice is transferred to the frost |
|---|
| 253 | is_h2o_perice(i) = .false. |
|---|
| 254 | h2o_frost4PCM(i,:) = h2o_frost4PCM(i,:) + h2o_ice(i,:) |
|---|
| 255 | h2o_ice(i,:) = 0._dp |
|---|
| 256 | end if |
|---|
| 257 | end do |
|---|
| 258 | |
|---|
| 259 | END SUBROUTINE build4PCM_perice |
|---|
| 260 | !======================================================================= |
|---|
| 261 | |
|---|
| 262 | !======================================================================= |
|---|
| 263 | SUBROUTINE evolve_co2ice(co2_ice,d_co2ice,zshift_surf) |
|---|
| 264 | !----------------------------------------------------------------------- |
|---|
| 265 | ! NAME |
|---|
| 266 | ! evolve_co2ice |
|---|
| 267 | ! |
|---|
| 268 | ! DESCRIPTION |
|---|
| 269 | ! Compute the evolution of CO2 ice over one year. |
|---|
| 270 | ! |
|---|
| 271 | ! AUTHORS & DATE |
|---|
| 272 | ! R. Vandemeulebrouck |
|---|
| 273 | ! L. Lange |
|---|
| 274 | ! JB Clement, 2023-2025 |
|---|
| 275 | ! |
|---|
| 276 | ! NOTES |
|---|
| 277 | ! |
|---|
| 278 | !----------------------------------------------------------------------- |
|---|
| 279 | |
|---|
| 280 | ! DEPENDENCIES |
|---|
| 281 | ! ------------ |
|---|
| 282 | use geometry, only: ngrid, nslope |
|---|
| 283 | use evolution, only: dt |
|---|
| 284 | use display, only: print_msg |
|---|
| 285 | |
|---|
| 286 | ! DECLARATION |
|---|
| 287 | ! ----------- |
|---|
| 288 | implicit none |
|---|
| 289 | |
|---|
| 290 | ! ARGUMENTS |
|---|
| 291 | ! --------- |
|---|
| 292 | real(dp), dimension(:,:), intent(inout) :: co2_ice |
|---|
| 293 | real(dp), dimension(:,:), intent(inout) :: d_co2ice |
|---|
| 294 | real(dp), dimension(:,:), intent(out) :: zshift_surf |
|---|
| 295 | |
|---|
| 296 | ! LOCAL VARIABLES |
|---|
| 297 | ! --------------- |
|---|
| 298 | real(dp), dimension(ngrid,nslope) :: co2_ice_old ! Old density of CO2 ice |
|---|
| 299 | |
|---|
| 300 | ! CODE |
|---|
| 301 | ! ---- |
|---|
| 302 | ! Evolution of CO2 ice for each physical point |
|---|
| 303 | call print_msg('> Evolution of CO2 ice') |
|---|
| 304 | |
|---|
| 305 | co2_ice_old = co2_ice |
|---|
| 306 | co2_ice = co2_ice + d_co2ice*dt |
|---|
| 307 | where (co2_ice < 0._dp) |
|---|
| 308 | co2_ice = 0._dp |
|---|
| 309 | d_co2ice = -co2_ice_old/dt |
|---|
| 310 | end where |
|---|
| 311 | |
|---|
| 312 | zshift_surf = d_co2ice*dt/rho_co2ice |
|---|
| 313 | |
|---|
| 314 | END SUBROUTINE evolve_co2ice |
|---|
| 315 | !======================================================================= |
|---|
| 316 | |
|---|
| 317 | !======================================================================= |
|---|
| 318 | SUBROUTINE evolve_h2oice(delta_h2o_ads,delta_icetable,h2o_ice,d_h2oice,zshift_surf,stopcrit) |
|---|
| 319 | !----------------------------------------------------------------------- |
|---|
| 320 | ! NAME |
|---|
| 321 | ! evolve_h2oice |
|---|
| 322 | ! |
|---|
| 323 | ! DESCRIPTION |
|---|
| 324 | ! Compute the evolution of H2O ice over one year. |
|---|
| 325 | ! |
|---|
| 326 | ! AUTHORS & DATE |
|---|
| 327 | ! R. Vandemeulebrouck |
|---|
| 328 | ! L. Lange |
|---|
| 329 | ! JB Clement, 2023-2025 |
|---|
| 330 | ! |
|---|
| 331 | ! NOTES |
|---|
| 332 | ! |
|---|
| 333 | !----------------------------------------------------------------------- |
|---|
| 334 | |
|---|
| 335 | ! DEPENDENCIES |
|---|
| 336 | ! ------------ |
|---|
| 337 | use geometry, only: ngrid, nslope |
|---|
| 338 | use evolution, only: dt |
|---|
| 339 | use stopping_crit, only: stopping_crit_h2o, stopFlags |
|---|
| 340 | use display, only: print_msg |
|---|
| 341 | |
|---|
| 342 | ! DECLARATION |
|---|
| 343 | ! ----------- |
|---|
| 344 | implicit none |
|---|
| 345 | |
|---|
| 346 | ! ARGUMENTS |
|---|
| 347 | ! --------- |
|---|
| 348 | real(dp), dimension(:), intent(in) :: delta_h2o_ads, delta_icetable |
|---|
| 349 | real(dp), dimension(:,:), intent(inout) :: h2o_ice, d_h2oice |
|---|
| 350 | type(stopFlags), intent(inout) :: stopcrit |
|---|
| 351 | real(dp), dimension(:,:), intent(out) :: zshift_surf |
|---|
| 352 | |
|---|
| 353 | ! LOCAL VARIABLES |
|---|
| 354 | ! --------------- |
|---|
| 355 | real(qp) :: S_atm_2_h2o, S_h2o_2_atm, S_atm_2_h2oice, S_h2oice_2_atm ! H2O balance variables |
|---|
| 356 | real(dp), dimension(ngrid,nslope) :: d_h2oice_new ! Tendencies computed to keep balance |
|---|
| 357 | |
|---|
| 358 | ! CODE |
|---|
| 359 | ! ---- |
|---|
| 360 | call print_msg('> Evolution of H2O ice') |
|---|
| 361 | |
|---|
| 362 | if (ngrid == 1) then ! In 1D |
|---|
| 363 | where (d_h2oice(:,:) < 0._dp .and. abs(d_h2oice(:,:)*dt) > h2o_ice(:,:)) ! Not enough ice to sublimate everything |
|---|
| 364 | ! We sublimate what we can |
|---|
| 365 | d_h2oice_new(:,:) = h2o_ice(:,:)/dt |
|---|
| 366 | ! It means the tendency is zero next time |
|---|
| 367 | d_h2oice(:,:) = 0._dp |
|---|
| 368 | else where |
|---|
| 369 | d_h2oice_new(:,:) = d_h2oice(:,:) |
|---|
| 370 | end where |
|---|
| 371 | else ! In 3D |
|---|
| 372 | call stopping_crit_h2o(delta_h2o_ads,delta_icetable,h2o_ice,d_h2oice,S_atm_2_h2o,S_h2o_2_atm,S_atm_2_h2oice,S_h2oice_2_atm,stopcrit) |
|---|
| 373 | if (stopcrit%stop_code() > 0) return |
|---|
| 374 | call balance_h2oice_reservoirs(S_atm_2_h2o,S_h2o_2_atm,S_atm_2_h2oice,S_h2oice_2_atm,h2o_ice,d_h2oice,d_h2oice_new) |
|---|
| 375 | end if |
|---|
| 376 | |
|---|
| 377 | h2o_ice = h2o_ice + d_h2oice_new*dt |
|---|
| 378 | zshift_surf = d_h2oice_new*dt/rho_h2oice |
|---|
| 379 | |
|---|
| 380 | END SUBROUTINE evolve_h2oice |
|---|
| 381 | !======================================================================= |
|---|
| 382 | |
|---|
| 383 | !======================================================================= |
|---|
| 384 | SUBROUTINE balance_h2oice_reservoirs(S_atm_2_h2o,S_h2o_2_atm,S_atm_2_h2oice,S_h2oice_2_atm,h2o_ice,d_h2oice,d_h2oice_new) |
|---|
| 385 | !----------------------------------------------------------------------- |
|---|
| 386 | ! NAME |
|---|
| 387 | ! balance_h2oice_reservoirs |
|---|
| 388 | ! |
|---|
| 389 | ! DESCRIPTION |
|---|
| 390 | ! Balance H2O flux from and into atmosphere across reservoirs. |
|---|
| 391 | ! |
|---|
| 392 | ! AUTHORS & DATE |
|---|
| 393 | ! JB Clement, 2025 |
|---|
| 394 | ! |
|---|
| 395 | ! NOTES |
|---|
| 396 | ! |
|---|
| 397 | !----------------------------------------------------------------------- |
|---|
| 398 | |
|---|
| 399 | ! DEPENDENCIES |
|---|
| 400 | ! ------------ |
|---|
| 401 | use evolution, only: dt |
|---|
| 402 | use geometry, only: ngrid, nslope |
|---|
| 403 | |
|---|
| 404 | ! DECLARATION |
|---|
| 405 | ! ----------- |
|---|
| 406 | implicit none |
|---|
| 407 | |
|---|
| 408 | ! ARGUMENTS |
|---|
| 409 | ! --------- |
|---|
| 410 | real(qp), intent(in) :: S_atm_2_h2o, S_h2o_2_atm, S_atm_2_h2oice, S_h2oice_2_atm |
|---|
| 411 | real(dp), dimension(:,:), intent(in) :: h2o_ice |
|---|
| 412 | real(dp), dimension(:,:), intent(inout) :: d_h2oice |
|---|
| 413 | real(dp), dimension(:,:), intent(out) :: d_h2oice_new |
|---|
| 414 | |
|---|
| 415 | ! LOCAL VARIABLES |
|---|
| 416 | ! --------------- |
|---|
| 417 | integer(di) :: i, islope |
|---|
| 418 | real(qp) :: S_target, S_target_subl_h2oice, S_target_cond_h2oice, S_ghostice ! Balance variables |
|---|
| 419 | real(dp) :: d_target |
|---|
| 420 | |
|---|
| 421 | ! CODE |
|---|
| 422 | ! ---- |
|---|
| 423 | S_target = (S_atm_2_h2o + S_h2o_2_atm)/2._dp |
|---|
| 424 | S_target_cond_h2oice = S_atm_2_h2oice + S_target - S_atm_2_h2o |
|---|
| 425 | S_target_subl_h2oice = S_h2oice_2_atm + S_target - S_h2o_2_atm |
|---|
| 426 | |
|---|
| 427 | d_h2oice_new = 0._dp |
|---|
| 428 | S_ghostice = 0._qp |
|---|
| 429 | do i = 1,ngrid |
|---|
| 430 | do islope = 1,nslope |
|---|
| 431 | if (d_h2oice(i,islope) > 0._dp) then ! Condensing |
|---|
| 432 | d_h2oice_new(i,islope) = d_h2oice(i,islope)*real(S_target_cond_h2oice/S_atm_2_h2oice,dp) |
|---|
| 433 | else if (d_h2oice(i,islope) < 0._dp) then ! Sublimating |
|---|
| 434 | d_target = d_h2oice(i,islope)*real(S_target_subl_h2oice/S_h2oice_2_atm,dp) |
|---|
| 435 | if (abs(d_target*dt) <= h2o_ice(i,islope)) then ! Enough ice to sublimate everything |
|---|
| 436 | d_h2oice_new(i,islope) = d_target |
|---|
| 437 | else ! Not enough ice to sublimate everything |
|---|
| 438 | ! We sublimate what we can |
|---|
| 439 | d_h2oice_new(i,islope) = h2o_ice(i,islope)/dt |
|---|
| 440 | ! It means the tendency is zero next time |
|---|
| 441 | d_h2oice(i,islope) = 0._dp |
|---|
| 442 | ! We compute the amount of H2O ice that we could not make sublimate |
|---|
| 443 | S_ghostice = S_ghostice + abs(d_target*dt) - h2o_ice(i,islope) |
|---|
| 444 | end if |
|---|
| 445 | end if |
|---|
| 446 | end do |
|---|
| 447 | end do |
|---|
| 448 | |
|---|
| 449 | ! We need to remove this ice unable to sublimate from places where ice condensed in order to keep balance |
|---|
| 450 | where (d_h2oice_new > 0._dp) d_h2oice_new = d_h2oice_new*real((S_target_cond_h2oice - S_ghostice)/S_target_cond_h2oice,dp) |
|---|
| 451 | |
|---|
| 452 | END SUBROUTINE balance_h2oice_reservoirs |
|---|
| 453 | !======================================================================= |
|---|
| 454 | |
|---|
| 455 | END MODULE surf_ice |
|---|