PISM, A Parallel Ice Sheet Model 2.3.0-79cae578d committed by Constantine Khrulev on 2026-03-22
Loading...
Searching...
No Matches
Routing.cc
Go to the documentation of this file.
1// Copyright (C) 2012-2023, 2025, 2026 PISM Authors
2//
3// This file is part of PISM.
4//
5// PISM is free software; you can redistribute it and/or modify it under the
6// terms of the GNU General Public License as published by the Free Software
7// Foundation; either version 3 of the License, or (at your option) any later
8// version.
9//
10// PISM is distributed in the hope that it will be useful, but WITHOUT ANY
11// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
12// FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
13// details.
14//
15// You should have received a copy of the GNU General Public License
16// along with PISM; if not, write to the Free Software
17// Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
18
19#include <cassert>
20#include <limits>
21
22#include "Hydrology.hh"
23#include "pism/hydrology/Routing.hh"
24#include "pism/util/array/CellType.hh"
25
26#include "pism/util/error_handling.hh"
27
28#include "pism/util/pism_utilities.hh"
29#include "pism/util/Vars.hh"
30#include "pism/geometry/Geometry.hh"
31#include "pism/util/Profiling.hh"
32#include "pism/util/Context.hh"
33#include "pism/util/Logger.hh"
34#include "pism/util/io/IO_Flags.hh"
35
36namespace pism {
37namespace hydrology {
38
39namespace diagnostics {
40
41//! \brief Reports the pressure of the transportable water in the subglacial layer.
42class BasalWaterPressure : public Diag<Routing>
43{
44public:
46 : Diag<Routing>(m) {
47 m_vars = { { m_sys, "bwp", *m_grid } };
48 m_vars[0].long_name("pressure of transportable water in subglacial layer").units("Pa");
49 }
50
51protected:
52 virtual std::shared_ptr<array::Array> compute_impl() const {
53 auto result = allocate<array::Scalar>("bwp");
54
55 result->copy_from(model->subglacial_water_pressure());
56
57 return result;
58 }
59};
60
61
62//! \brief Reports the pressure of the transportable water in the subglacial layer as a
63//! fraction of the overburden pressure.
64class RelativeBasalWaterPressure : public Diag<Routing>
65{
66public:
68 : Diag<Routing>(m) {
69 m_vars = { { m_sys, "bwprel", *m_grid } };
70 m_vars[0]
71 .long_name(
72 "pressure of transportable water in subglacial layer as fraction of the overburden pressure")
73 .units("1");
74 m_vars[0]["_FillValue"] = {m_fill_value};
75 }
76
77protected:
78 virtual std::shared_ptr<array::Array> compute_impl() const {
79 double fill_value = m_fill_value;
80
81 auto result = allocate<array::Scalar>("bwprel");
82
83 const array::Scalar
84 &P = model->subglacial_water_pressure(),
85 &Po = model->overburden_pressure();
86
87 array::AccessScope list{result.get(), &Po, &P};
88 for (auto p : m_grid->points()) {
89 const int i = p.i(), j = p.j();
90
91 if (Po(i,j) > 0.0) {
92 (*result)(i,j) = P(i, j) / Po(i,j);
93 } else {
94 (*result)(i,j) = fill_value;
95 }
96 }
97
98 return result;
99 }
100};
101
102
103//! \brief Reports the effective pressure of the transportable water in the subglacial
104//! layer, that is, the overburden pressure minus the pressure.
105class EffectiveBasalWaterPressure : public Diag<Routing>
106{
107public:
109 : Diag<Routing>(m) {
110 m_vars = { { m_sys, "effbwp", *m_grid } };
111 m_vars[0]
112 .long_name("effective pressure of transportable water in subglacial layer"
113 " (overburden pressure minus water pressure)")
114 .units("Pa");
115 }
116
117protected:
118 virtual std::shared_ptr<array::Array> compute_impl() const {
119
120 auto result = allocate<array::Scalar>("effbwp");
121
122 const array::Scalar
123 &P = model->subglacial_water_pressure(),
124 &Po = model->overburden_pressure();
125
126 array::AccessScope list{&Po, &P, result.get()};
127
128 for (auto p : m_grid->points()) {
129 const int i = p.i(), j = p.j();
130
131 (*result)(i, j) = Po(i, j) - P(i, j);
132 }
133
134 return result;
135 }
136};
137
138
139//! \brief Report the wall melt rate from dissipation of the potential energy of the
140//! transportable water.
141class WallMelt : public Diag<Routing>
142{
143public:
145 : Diag<Routing>(m) {
146 m_vars = { { m_sys, "wallmelt", *m_grid } };
147 m_vars[0]
148 .long_name("wall melt into subglacial hydrology layer from (turbulent)"
149 " dissipation of energy in transportable water")
150 .units("m s^-1")
151 .output_units("m year^-1");
152 }
153
154protected:
155 virtual std::shared_ptr<array::Array> compute_impl() const {
156 auto result = allocate<array::Scalar>("wallmelt");
157
158 const array::Scalar &bed_elevation = *m_grid->variables().get_2d_scalar("bedrock_altitude");
159
160 wall_melt(*model, bed_elevation, *result);
161 return result;
162 }
163};
164
165//! @brief Diagnostically reports the staggered-grid components of the velocity of the
166//! water in the subglacial layer.
167class BasalWaterVelocity : public Diag<Routing>
168{
169public:
171 : Diag<Routing>(m) {
172 m_vars = { { m_sys, "bwatvel[0]", *m_grid }, { m_sys, "bwatvel[1]", *m_grid } };
173 m_vars[0].long_name("velocity of water in subglacial layer, i-offset").units("m s^-1");
174 m_vars[1].long_name("velocity of water in subglacial layer, j-offset").units("m s^-1");
175 }
176protected:
177 virtual std::shared_ptr<array::Array> compute_impl() const {
178 auto result = allocate<array::Staggered>("bwatvel");
179
180 result->copy_from(model->velocity_staggered());
181
182 return result;
183 }
184};
185
186//! Compute the hydraulic potential.
187/*!
188 Computes \f$\psi = P + \rho_w g (b + W)\f$.
189*/
191 const array::Scalar &P,
192 const array::Scalar &sea_level,
193 const array::Scalar &bed,
194 const array::Scalar &ice_thickness,
195 array::Scalar &result) {
196
197 auto grid = result.grid();
198
199 auto config = grid->ctx()->config();
200
201 double
202 ice_density = config->get_number("constants.ice.density"),
203 sea_water_density = config->get_number("constants.sea_water.density"),
204 C = ice_density / sea_water_density,
205 rg = (config->get_number("constants.fresh_water.density") *
206 config->get_number("constants.standard_gravity"));
207
208 array::AccessScope list{&P, &W, &sea_level, &ice_thickness, &bed, &result};
209
210 for (auto p : grid->points()) {
211 const int i = p.i(), j = p.j();
212
213 double b = std::max(bed(i, j), sea_level(i, j) - C * ice_thickness(i, j));
214
215 result(i, j) = P(i, j) + rg * (b + W(i, j));
216 }
217}
218
219/*! @brief Report hydraulic potential in the subglacial hydrology system */
220class HydraulicPotential : public Diag<Routing>
221{
222public:
224 m_vars = { { m_sys, "hydraulic_potential", *m_grid } };
225 m_vars[0].long_name("hydraulic potential in the subglacial hydrology system").units("Pa");
226 }
227
228protected:
229 std::shared_ptr<array::Array> compute_impl() const {
230
231 auto result = allocate<array::Scalar>("hydraulic_potential");
232
233 const auto &sea_level = *m_grid->variables().get_2d_scalar("sea_level");
234 const auto &bed_elevation = *m_grid->variables().get_2d_scalar("bedrock_altitude");
235 const auto &ice_thickness = *m_grid->variables().get_2d_scalar("land_ice_thickness");
236
237 hydraulic_potential(model->subglacial_water_thickness(),
238 model->subglacial_water_pressure(),
239 sea_level,
240 bed_elevation,
241 ice_thickness,
242 *result);
243
244 return result;
245 }
246};
247
248} // end of namespace diagnostics
249
250Routing::Routing(std::shared_ptr<const Grid> grid)
251 : Hydrology(grid),
252 m_Qstag(grid, "advection_flux"),
253 m_Qstag_average(grid, "cumulative_advection_flux"),
254 m_Vstag(grid, "water_velocity"),
255 m_Wstag(grid, "W_staggered"),
256 m_Kstag(grid, "K_staggered"),
257 m_Wnew(grid, "W_new"),
258 m_Wtillnew(grid, "Wtill_new"),
259 m_R(grid, "potential_workspace"), /* box stencil used */
260 m_dx(grid->dx()),
261 m_dy(grid->dy()),
262 m_bottom_surface(grid, "ice_bottom_surface_elevation") {
263
264 m_rg = (m_config->get_number("constants.fresh_water.density") *
265 m_config->get_number("constants.standard_gravity"));
266
268 .long_name("cell face-centered (staggered) components of advective subglacial water flux")
269 .units("m^2 s^-1");
270
272 .long_name("average (over time) advection flux on the staggered grid")
273 .units("m^2 s^-1");
274
276 .long_name(
277 "cell face-centered (staggered) components of water velocity in subglacial water layer")
278 .units("m s^-1");
279
280 // auxiliary variables which NEED ghosts
282 .long_name("cell face-centered (staggered) values of water layer thickness")
283 .units("m");
284 m_Wstag.metadata()["valid_min"] = { 0.0 };
285
287 .long_name("cell face-centered (staggered) values of nonlinear conductivity");
288 m_Kstag.metadata()["valid_min"] = { 0.0 };
289
290 m_R.metadata(0)
291 .long_name("work space for modeled subglacial water hydraulic potential")
292 .units("Pa");
293
294 // temporaries during update; do not need ghosts
296 .long_name("new thickness of transportable subglacial water layer during update")
297 .units("m");
298 m_Wnew.metadata()["valid_min"] = { 0.0 };
299
301 .long_name("new thickness of till (subglacial) water layer during update")
302 .units("m");
303 m_Wtillnew.metadata()["valid_min"] = { 0.0 };
304
305 {
306 double alpha = m_config->get_number("hydrology.thickness_power_in_flux");
307 if (alpha < 1.0) {
309 "alpha = %f < 1 which is not allowed", alpha);
310 }
311
312 if (m_config->get_number("hydrology.tillwat_max") < 0.0) {
314 "hydrology::Routing: hydrology.tillwat_max is negative.\n"
315 "This is not allowed.");
316 }
317 }
318}
319
321 m_log->message(2,
322 "* Initializing the routing subglacial hydrology model ...\n");
323
324 if (m_config->get_flag("hydrology.routing.include_floating_ice")) {
325 m_log->message(2, " ... routing subglacial water under grounded and floating ice.\n");
326 } else {
327 m_log->message(2, " ... routing subglacial water under grounded ice only.\n");
328 }
329}
330
331void Routing::restart_impl(const File &input_file, int record) {
332 Hydrology::restart_impl(input_file, record);
333
334 m_W.read(input_file, record);
335
336 regrid("Hydrology", m_W);
337}
338
339void Routing::bootstrap_impl(const File &input_file,
340 const array::Scalar &ice_thickness) {
341 Hydrology::bootstrap_impl(input_file, ice_thickness);
342
343 double bwat_default = m_config->get_number("bootstrapping.defaults.bwat");
344 m_W.regrid(input_file, io::Default(bwat_default));
345
346 regrid("Hydrology", m_W);
347}
348
350 const array::Scalar &W,
351 const array::Scalar &P) {
352 Hydrology::init_impl(W_till, W, P);
353
354 m_W.copy_from(W);
355}
356
357std::set<VariableMetadata> Routing::state_impl() const {
359}
360
361void Routing::write_state_impl(const OutputFile &output) const {
363 m_W.write(output);
364}
365
366//! Returns the (trivial) overburden pressure as the pressure of the transportable water,
367//! because this is the model.
371
373 return m_Vstag;
374}
375
376
377//! Average the regular grid water thickness to values at the center of cell edges.
378/*! Uses mask values to avoid averaging using water thickness values from
379 either ice-free or floating areas. */
381 const array::CellType1 &mask,
382 array::Staggered &result) {
383
384 bool include_floating = m_config->get_flag("hydrology.routing.include_floating_ice");
385
386 array::AccessScope list{ &mask, &W, &result };
387
388 for (auto p : m_grid->points()) {
389 const int i = p.i(), j = p.j();
390
391 if (include_floating) {
392 // east
393 if (mask.icy(i, j)) {
394 result(i, j, 0) = mask.icy(i + 1, j) ? 0.5 * (W(i, j) + W(i + 1, j)) : W(i, j);
395 } else {
396 result(i, j, 0) = mask.icy(i + 1, j) ? W(i + 1, j) : 0.0;
397 }
398 // north
399 if (mask.icy(i, j)) {
400 result(i, j, 1) = mask.icy(i, j + 1) ? 0.5 * (W(i, j) + W(i, j + 1)) : W(i, j);
401 } else {
402 result(i, j, 1) = mask.icy(i, j + 1) ? W(i, j + 1) : 0.0;
403 }
404 } else {
405 // east
406 if (mask.grounded_ice(i, j)) {
407 result(i, j, 0) = mask.grounded_ice(i + 1, j) ? 0.5 * (W(i, j) + W(i + 1, j)) : W(i, j);
408 } else {
409 result(i, j, 0) = mask.grounded_ice(i + 1, j) ? W(i + 1, j) : 0.0;
410 }
411 // north
412 if (mask.grounded_ice(i, j)) {
413 result(i, j, 1) = mask.grounded_ice(i, j + 1) ? 0.5 * (W(i, j) + W(i, j + 1)) : W(i, j);
414 } else {
415 result(i, j, 1) = mask.grounded_ice(i, j + 1) ? W(i, j + 1) : 0.0;
416 }
417 }
418 }
419
420 result.update_ghosts();
421}
422
423
424//! Compute the nonlinear conductivity at the center of cell edges.
425/*!
426 Computes
427
428 \f[ K = K(W, \nabla P, \nabla b) = k W^{\alpha-1} |\nabla R|^{\beta-2} \f]
429
430 on the staggered grid, where \f$R = P+\rho_w g b\f$. We denote
431
432 \f[ \Pi = |\nabla R|^2 \f]
433
434 internally; this is computed on a staggered grid by a Mahaffy-like ([@ref Mahaffy])
435 scheme. This requires \f$R\f$ to be defined on a box stencil of width 1.
436
437 Also returns the maximum over all staggered points of \f$ K W \f$.
438*/
440 const array::Scalar &P,
441 const array::Scalar &bed_elevation,
442 array::Staggered &result,
443 double &KW_max) const {
444 const double
445 k = m_config->get_number("hydrology.hydraulic_conductivity"),
446 alpha = m_config->get_number("hydrology.thickness_power_in_flux"),
447 beta = m_config->get_number("hydrology.gradient_power_in_flux"),
448 betapow = (beta - 2.0) / 2.0;
449
450 array::AccessScope list({&result, &W});
451
452 KW_max = 0.0;
453
454 if (beta != 2.0) {
455 // Put the squared norm of the gradient of the simplified hydrolic potential (Pi) in
456 // "result"
457 //
458 // FIXME: we don't need to re-compute this during every hydrology time step: the
459 // simplified hydrolic potential does not depend on the water amount and can be
460 // computed *once* in update_impl(), before entering the time-stepping loop
461 {
462 // R <-- P + rhow g b
463 P.add(m_rg, bed_elevation, m_R); // yes, it updates ghosts
464
465 list.add(m_R);
466 for (auto p : m_grid->points()) {
467 const int i = p.i(), j = p.j();
468
469 double dRdx, dRdy;
470 dRdx = (m_R(i + 1, j) - m_R(i, j)) / m_dx;
471 dRdy = (m_R(i + 1, j + 1) + m_R(i, j + 1) - m_R(i + 1, j - 1) - m_R(i, j - 1)) / (4.0 * m_dy);
472 result(i, j, 0) = dRdx * dRdx + dRdy * dRdy;
473
474 dRdx = (m_R(i + 1, j + 1) + m_R(i + 1, j) - m_R(i - 1, j + 1) - m_R(i - 1, j)) / (4.0 * m_dx);
475 dRdy = (m_R(i, j + 1) - m_R(i, j)) / m_dy;
476 result(i, j, 1) = dRdx * dRdx + dRdy * dRdy;
477 }
478 }
479
480 // We regularize negative power |\grad psi|^{beta-2} by adding eps because large
481 // head gradient might be 10^7 Pa per 10^4 m or 10^3 Pa/m.
482 const double eps = beta < 2.0 ? 1.0 : 0.0;
483
484 for (auto p : m_grid->points()) {
485 const int i = p.i(), j = p.j();
486
487 for (int o = 0; o < 2; ++o) {
488 const double Pi = result(i, j, o);
489
490 // FIXME: same as Pi above: we don't need to re-compute this each time we make a
491 // short hydrology time step
492 double B = pow(Pi + eps * eps, betapow);
493
494 result(i, j, o) = k * pow(W(i, j, o), alpha - 1.0) * B;
495
496 KW_max = std::max(KW_max, result(i, j, o) * W(i, j, o));
497 }
498 }
499 } else {
500 for (auto p : m_grid->points()) {
501 const int i = p.i(), j = p.j();
502
503 for (int o = 0; o < 2; ++o) {
504 result(i, j, o) = k * pow(W(i, j, o), alpha - 1.0);
505
506 KW_max = std::max(KW_max, result(i, j, o) * W(i, j, o));
507 }
508 }
509 }
510
511 KW_max = GlobalMax(m_grid->com, KW_max);
512
513 result.update_ghosts();
514}
515
516
517//! Compute the wall melt rate which comes from (turbulent) dissipation of flow energy.
518/*!
519 This code fills `result` with
520 \f[ \frac{m_{wall}}{\rho_w} = - \frac{1}{L \rho_w} \mathbf{q} \cdot \nabla \psi = \left(\frac{k}{L \rho_w}\right) W^\alpha |\nabla R|^\beta \f]
521 where \f$R = P+\rho_w g b\f$.
522
523 Note that conductivity_staggered() computes the related quantity
524 \f$K = k W^{\alpha-1} |\nabla R|^{\beta-2}\f$ on the staggered grid, but
525 contriving to reuse that code would be inefficient because of the
526 staggered-versus-regular change.
527
528 At the current state of the code, this is a diagnostic calculation only.
529*/
530void wall_melt(const Routing &model,
531 const array::Scalar &bed_elevation,
532 array::Scalar &result) {
533
534 auto grid = result.grid();
535
536 auto config = grid->ctx()->config();
537
538 const double
539 k = config->get_number("hydrology.hydraulic_conductivity"),
540 L = config->get_number("constants.fresh_water.latent_heat_of_fusion"),
541 alpha = config->get_number("hydrology.thickness_power_in_flux"),
542 beta = config->get_number("hydrology.gradient_power_in_flux"),
543 g = config->get_number("constants.standard_gravity"),
544 rhow = config->get_number("constants.fresh_water.density"),
545 rg = rhow * g,
546 CC = k / (L * rhow);
547
548 // FIXME: could be scaled with overall factor hydrology_coefficient_wall_melt ?
549 if (alpha < 1.0) {
551 "alpha = %f < 1 which is not allowed", alpha);
552 }
553
554 array::Scalar1 R(grid, "R");
555
556 // R <-- P + rhow g b
557 model.subglacial_water_pressure().add(rg, bed_elevation, R);
558 // yes, it updates ghosts
559
560 array::Scalar1 W(grid, "W");
562
563 array::AccessScope list{&R, &W, &result};
564
565 double dx = grid->dx();
566 double dy = grid->dy();
567
568 for (auto p : grid->points()) {
569 const int i = p.i(), j = p.j();
570 double dRdx, dRdy;
571
572 if (W(i, j) > 0.0) {
573 dRdx = 0.0;
574 if (W(i + 1, j) > 0.0) {
575 dRdx = (R(i + 1, j) - R(i, j)) / (2.0 * dx);
576 }
577 if (W(i - 1, j) > 0.0) {
578 dRdx += (R(i, j) - R(i - 1, j)) / (2.0 * dx);
579 }
580 dRdy = 0.0;
581 if (W(i, j + 1) > 0.0) {
582 dRdy = (R(i, j + 1) - R(i, j)) / (2.0 * dy);
583 }
584 if (W(i, j - 1) > 0.0) {
585 dRdy += (R(i, j) - R(i, j - 1)) / (2.0 * dy);
586 }
587 result(i, j) = CC * pow(W(i, j), alpha) * pow(dRdx * dRdx + dRdy * dRdy, beta/2.0);
588 } else {
589 result(i, j) = 0.0;
590 }
591 }
592}
593
594
595//! Get the advection velocity V at the center of cell edges.
596/*!
597 Computes the advection velocity @f$\mathbf{V}@f$ on the staggered
598 (edge-centered) grid. If V = (u, v) in components then we have
599 <code> result(i, j, 0) = u(i+1/2, j) </code> and
600 <code> result(i, j, 1) = v(i, j+1/2) </code>
601
602 The advection velocity is given by the formula
603
604 @f[ \mathbf{V} = - K \left(\nabla P + \rho_w g \nabla b\right) @f]
605
606 where @f$\mathbf{V}@f$ is the water velocity, @f$P@f$ is the water
607 pressure, and @f$b@f$ is the bedrock elevation.
608
609 If the corresponding staggered grid value of the water thickness is zero then that
610 component of V is set to zero. This does not change the flux value (which would be zero
611 anyway) but it does provide the correct max velocity in the CFL calculation. We assume
612 bed has valid ghosts.
613*/
615 const array::Scalar &pressure,
616 const array::Scalar &bed,
617 const array::Staggered &K,
618 const array::Scalar1 *no_model_mask,
619 array::Staggered &result) const {
620 array::Scalar &P = m_R;
621 P.copy_from(pressure); // yes, it updates ghosts
622
623 array::AccessScope list{&P, &W, &K, &bed, &result};
624
625 for (auto p : m_grid->points()) {
626 const int i = p.i(), j = p.j();
627
628 if (W(i, j, 0) > 0.0) {
629 double
630 P_x = (P(i + 1, j) - P(i, j)) / m_dx,
631 b_x = (bed(i + 1, j) - bed(i, j)) / m_dx;
632 result(i, j, 0) = - K(i, j, 0) * (P_x + m_rg * b_x);
633 } else {
634 result(i, j, 0) = 0.0;
635 }
636
637 if (W(i, j, 1) > 0.0) {
638 double
639 P_y = (P(i, j + 1) - P(i, j)) / m_dy,
640 b_y = (bed(i, j + 1) - bed(i, j)) / m_dy;
641 result(i, j, 1) = - K(i, j, 1) * (P_y + m_rg * b_y);
642 } else {
643 result(i, j, 1) = 0.0;
644 }
645 }
646
647 if (no_model_mask) {
648 list.add(*no_model_mask);
649
650 for (auto p : m_grid->points()) {
651 const int i = p.i(), j = p.j();
652
653 auto M = no_model_mask->star(i, j);
654
655 if (M.c or M.e) {
656 result(i, j, 0) = 0.0;
657 }
658
659 if (M.c or M.n) {
660 result(i, j, 1) = 0.0;
661 }
662 }
663 }
664}
665
666
667//! Compute Q = V W at edge-centers (staggered grid) by first-order upwinding.
668/*!
669 The field W must have valid ghost values, but V does not need them.
670
671 FIXME: This could be re-implemented using the Koren (1993) flux-limiter.
672*/
674 const array::Scalar &W,
675 array::Staggered &result) const {
676 array::AccessScope list{&W, &V, &result};
677
678 assert(W.stencil_width() >= 1);
679
680 for (auto p : m_grid->points()) {
681 const int i = p.i(), j = p.j();
682
683 result(i, j, 0) = V(i, j, 0) * (V(i, j, 0) >= 0.0 ? W(i, j) : W(i + 1, j));
684 result(i, j, 1) = V(i, j, 1) * (V(i, j, 1) >= 0.0 ? W(i, j) : W(i, j + 1));
685 }
686
687 result.update_ghosts();
688}
689
690/*!
691 * See equation (51) in Bueler and van Pelt.
692 */
693double Routing::max_timestep_W_diff(double KW_max) const {
694 if (KW_max > 0.0) {
695 double D_max = m_rg * KW_max;
696 double result = 1.0 / (m_dx * m_dx) + 1.0 / (m_dy * m_dy);
697 return 0.25 / (D_max * result);
698 }
699 return std::numeric_limits<double>::infinity();
700}
701
702/*!
703 * See equation (50) in Bueler and van Pelt.
704 */
706 // V could be zero if P is constant and bed is flat
707 auto tmp = absmax(m_Vstag);
708
709 // add a safety margin
710 double alpha = 0.95;
711 double eps = 1e-6;
712
713 return alpha * 0.5 / (tmp[0]/m_dx + tmp[1]/m_dy + eps);
714}
715
716
717//! The computation of Wtillnew, called by update().
718/*!
719 Does a step of the trivial integration
720 \f[ \frac{\partial W_{till}}{\partial t} = \frac{m}{\rho_w} - C\f]
721
722 where \f$C=\f$`hydrology_tillwat_decay_rate`. Enforces bounds
723 \f$0 \le W_{till} \le W_{till}^{max}\f$ where the upper bound is
724 `hydrology_tillwat_max`. Here \f$m/\rho_w\f$ is `total_input`.
725
726 Compare hydrology::NullTransport::update_impl().
727
728 The current code is not quite "code duplication" because the code here:
729
730 1. computes `Wtill_new` instead of updating `Wtill` in place;
731 2. uses time steps determined by the rest of the hydrology::Routing model;
732 3. does not check mask because the enforce_bounds() call addresses that.
733
734 Otherwise this is the same physical model with the same configurable parameters.
735*/
737 const array::Scalar &Wtill,
738 const array::Scalar &surface_input_rate,
739 const array::Scalar &basal_melt_rate,
740 array::Scalar &Wtill_new) {
741 const double
742 tillwat_max = m_config->get_number("hydrology.tillwat_max"),
743 C = m_config->get_number("hydrology.tillwat_decay_rate", "m / second");
744
745 array::AccessScope list{&Wtill, &Wtill_new, &basal_melt_rate};
746
747 bool add_surface_input = m_config->get_flag("hydrology.add_water_input_to_till_storage");
748 if (add_surface_input) {
750 }
751
752 for (auto p : m_grid->points()) {
753 const int i = p.i(), j = p.j();
754
755 double input_rate = basal_melt_rate(i, j);
756 if (add_surface_input) {
757 input_rate += surface_input_rate(i, j);
758 }
759
760 Wtill_new(i, j) = clip(Wtill(i, j) + dt * (input_rate - C),
761 0.0, tillwat_max);
762 }
763}
764
766 const array::Scalar1 &W,
767 const array::Staggered1 &Wstag,
768 const array::Staggered1 &K,
769 const array::Staggered1 &Q,
770 array::Scalar &result) {
771 const double
772 wux = 1.0 / (m_dx * m_dx),
773 wuy = 1.0 / (m_dy * m_dy);
774
775 array::AccessScope list{&W, &Wstag, &K, &Q, &result};
776
777 for (auto p : m_grid->points()) {
778 const int i = p.i(), j = p.j();
779
780 auto q = Q.star(i, j);
781 const double divQ = (q.e - q.w) / m_dx + (q.n - q.s) / m_dy;
782
783 auto k = K.star(i, j);
784 auto ws = Wstag.star(i, j);
785
786 const double
787 De = m_rg * k.e * ws.e,
788 Dw = m_rg * k.w * ws.w,
789 Dn = m_rg * k.n * ws.n,
790 Ds = m_rg * k.s * ws.s;
791
792 auto w = W.star(i, j);
793 const double diffW = (wux * (De * (w.e - w.c) - Dw * (w.c - w.w)) +
794 wuy * (Dn * (w.n - w.c) - Ds * (w.c - w.s)));
795
796 result(i, j) = dt * (- divQ + diffW);
797 }
798}
799
800
801//! The computation of Wnew, called by update().
802void Routing::update_W(double dt,
803 const array::Scalar &surface_input_rate,
804 const array::Scalar &basal_melt_rate,
805 const array::Scalar1 &W,
806 const array::Staggered1 &Wstag,
807 const array::Scalar &Wtill,
808 const array::Scalar &Wtill_new,
809 const array::Staggered1 &K,
810 const array::Staggered1 &Q,
811 array::Scalar &W_new) {
812
814
815 array::AccessScope list{&W, &Wtill, &Wtill_new, &surface_input_rate,
816 &basal_melt_rate, &m_flow_change_incremental, &W_new};
817
818 for (auto p : m_grid->points()) {
819 const int i = p.i(), j = p.j();
820
821 double input_rate = surface_input_rate(i, j) + basal_melt_rate(i, j);
822
823 double Wtill_change = Wtill_new(i, j) - Wtill(i, j);
824 W_new(i, j) = (W(i, j) + (dt * input_rate - Wtill_change) + m_flow_change_incremental(i, j));
825 }
826
829 m_input_change.add(dt, basal_melt_rate);
830}
831
832//! Update the model state variables W and Wtill by applying the subglacial hydrology model equations.
833/*!
834 Runs the hydrology model from time t to time t + dt. Here [t, dt]
835 is generally on the order of months to years. This hydrology model will take its
836 own shorter time steps, perhaps hours to weeks.
837
838 To update W = `bwat` we call update_W(), and to update Wtill = `tillwat` we
839 call update_Wtill().
840*/
841void Routing::update_impl(double t, double dt, const Inputs& inputs) {
842
844
845 double
846 ht = t,
847 hdt = 0.0;
848
849 const double
850 t_final = t + dt,
851 dt_max = m_config->get_number("hydrology.maximum_time_step", "seconds"),
852 tillwat_max = m_config->get_number("hydrology.tillwat_max");
853
854 m_Qstag_average.set(0.0);
855
856 // make sure W has valid ghosts before starting hydrology steps
858
859 unsigned int step_counter = 0;
860 for (; ht < t_final; ht += hdt) {
861 step_counter++;
862
863#if (Pism_DEBUG==1)
864 double huge_number = 1e6;
865 check_bounds(m_W, huge_number);
866 check_bounds(m_Wtill, tillwat_max);
867#endif
868
869 // updates ghosts of m_Wstag
871 inputs.geometry->cell_type,
872 m_Wstag);
873
874 double maxKW = 0.0;
875 // updates ghosts of m_Kstag
876 profiling().begin("routing_conductivity");
880 m_Kstag, maxKW);
881 profiling().end("routing_conductivity");
882
883 // ghosts of m_Vstag are not updated
884 profiling().begin("routing_velocity");
888 m_Kstag,
889 inputs.no_model_mask,
890 m_Vstag);
891 profiling().end("routing_velocity");
892
893 // to get Q, W needs valid ghosts (ghosts of m_Vstag are not used)
894 // updates ghosts of m_Qstag
895 profiling().begin("routing_flux");
897 profiling().end("routing_flux");
898
900
901 {
902 const double
903 dt_cfl = max_timestep_W_cfl(),
904 dt_diff_w = max_timestep_W_diff(maxKW);
905
906 hdt = std::min(t_final - ht, dt_max);
907 hdt = std::min(hdt, dt_cfl);
908 hdt = std::min(hdt, dt_diff_w);
909 }
910
911 m_log->message(3, " hydrology step %05d, dt = %f s\n", step_counter, hdt);
912
913 // update Wtillnew from Wtill and input_rate
914 {
915 profiling().begin("routing_Wtill");
916 update_Wtill(hdt,
917 m_Wtill,
920 m_Wtillnew);
921 // remove water in ice-free areas and account for changes
923 inputs.no_model_mask,
924 0.0, // do not limit maximum thickness
925 tillwat_max, // till water thickness under the ocean
931 profiling().end("routing_Wtill");
932 }
933
934 // update Wnew from W, Wtill, Wtillnew, Wstag, Q, input_rate
935 // uses ghosts of m_W, m_Wstag, m_Qstag, m_Kstag
936 {
937 profiling().begin("routing_W");
938 update_W(hdt,
941 m_W, m_Wstag,
944 m_Wnew);
945 // remove water in ice-free areas and account for changes
947 inputs.no_model_mask,
948 0.0, // do not limit maximum thickness
949 0.0, // transportable water thickness under the ocean
950 m_Wnew,
955
956 // transfer new into old (updates ghosts of m_W)
958 profiling().end("routing_W");
959 }
960
961 // m_Wtill has no ghosts
963 } // end of the time-stepping loop
964
965 staggered_to_regular(inputs.geometry->cell_type, m_Qstag_average,
966 m_config->get_flag("hydrology.routing.include_floating_ice"),
967 m_Q);
968 m_Q.scale(1.0 / dt);
969
970 m_log->message(2,
971 " took %d hydrology sub-steps with average dt = %.6f years (%.3f s or %.3f hours)\n",
972 step_counter,
973 units::convert(m_sys, dt / step_counter, "seconds", "years"),
974 dt / step_counter,
975 (dt / step_counter) / 3600.0);
976}
977
978std::map<std::string, Diagnostic::Ptr> Routing::spatial_diagnostics_impl() const {
979 using namespace diagnostics;
980
981 DiagnosticList result = {
982 {"bwatvel", Diagnostic::Ptr(new BasalWaterVelocity(this))},
983 {"bwp", Diagnostic::Ptr(new BasalWaterPressure(this))},
984 {"bwprel", Diagnostic::Ptr(new RelativeBasalWaterPressure(this))},
985 {"effbwp", Diagnostic::Ptr(new EffectiveBasalWaterPressure(this))},
986 {"wallmelt", Diagnostic::Ptr(new WallMelt(this))},
987 {"hydraulic_potential", Diagnostic::Ptr(new HydraulicPotential(this))},
988 };
990}
991
992std::map<std::string, TSDiagnostic::Ptr> Routing::scalar_diagnostics_impl() const {
993 std::map<std::string, TSDiagnostic::Ptr> result = {
994 // FIXME: add mass-conservation diagnostics
995 };
996 return result;
997}
998
999} // end of namespace hydrology
1000} // end of namespace pism
const units::System::Ptr m_sys
unit system used by this component
Definition Component.hh:162
std::shared_ptr< const Config > m_config
configuration database used by this component
Definition Component.hh:160
const std::shared_ptr< const Grid > m_grid
grid used by this component
Definition Component.hh:158
void regrid(const std::string &module_name, array::Array &variable, RegriddingFlag flag=NO_REGRID_WITHOUT_REGRID_VARS)
Definition Component.cc:152
const Profiling & profiling() const
Definition Component.cc:115
std::shared_ptr< const Logger > m_log
logger (for easy access)
Definition Component.hh:164
const Routing * model
A template derived from Diagnostic, adding a "Model".
std::vector< VariableMetadata > m_vars
metadata corresponding to NetCDF variables
double m_fill_value
fill value (used often enough to justify storing it)
const units::System::Ptr m_sys
the unit system
std::shared_ptr< Diagnostic > Ptr
Definition Diagnostic.hh:67
std::shared_ptr< const Grid > m_grid
the grid
High-level PISM I/O class.
Definition File.hh:57
array::CellType2 cell_type
Definition Geometry.hh:55
void begin(const char *name) const
Definition Profiling.cc:75
void end(const char *name) const
Definition Profiling.cc:91
static RuntimeError formatted(const ErrorLocation &location, const char format[],...) __attribute__((format(printf
build a RuntimeError with a formatted message
VariableMetadata & long_name(const std::string &input)
VariableMetadata & units(const std::string &input)
double get_number(const std::string &name) const
Get a single-valued scalar attribute.
Makes sure that we call begin_access() and end_access() for all accessed array::Arrays.
Definition Array.hh:66
void copy_from(const Array2D< T > &source)
Definition Array2D.hh:101
stencils::Star< T > star(int i, int j) const
Definition Array2D.hh:107
void add(double alpha, const Array2D< T > &x)
Definition Array2D.hh:93
void read(const std::string &filename, unsigned int time)
Definition Array.cc:753
void write(const OutputFile &file) const
Definition Array.cc:859
void scale(double alpha)
Result: v <- v * alpha. Calls VecScale.
Definition Array.cc:227
std::shared_ptr< const Grid > grid() const
Definition Array.cc:134
void set(double c)
Result: v[j] <- c for all j.
Definition Array.cc:659
void regrid(const std::string &filename, io::Default default_value)
Definition Array.cc:758
void add(double alpha, const Array &x)
Result: v <- v + alpha * x. Calls VecAXPY.
Definition Array.cc:209
void update_ghosts()
Updates ghost points.
Definition Array.cc:645
VariableMetadata & metadata(unsigned int N=0)
Returns a reference to the VariableMetadata object containing metadata for the compoment N.
Definition Array.cc:467
unsigned int stencil_width() const
Get the stencil width of the current Array. Returns 0 if ghosts are not available.
Definition Array.cc:305
bool grounded_ice(int i, int j) const
Definition CellType.hh:46
bool icy(int i, int j) const
Definition CellType.hh:42
stencils::Star< double > star(int i, int j) const
Returns the values at interfaces of the cell i,j using the staggered grid.
Definition Staggered.hh:75
A class for storing and accessing internal staggered-grid 2D fields. Uses dof=2 storage....
Definition Staggered.hh:37
array::Scalar m_flow_change
Definition Hydrology.hh:196
array::Scalar m_flow_change_incremental
Definition Hydrology.hh:185
virtual void bootstrap_impl(const File &input_file, const array::Scalar &ice_thickness)
Definition Hydrology.cc:378
virtual std::map< std::string, Diagnostic::Ptr > spatial_diagnostics_impl() const
Definition Hydrology.cc:444
array::Scalar m_surface_input_rate
Definition Hydrology.hh:179
array::Scalar m_grounding_line_change
Definition Hydrology.hh:192
virtual std::set< VariableMetadata > state_impl() const
Definition Hydrology.cc:468
const array::Scalar & surface_input_rate() const
Definition Hydrology.cc:519
void enforce_bounds(const array::CellType &cell_type, const array::Scalar *no_model_mask, double max_thickness, double ocean_water_thickness, array::Scalar &water_thickness, array::Scalar &grounded_margin_change, array::Scalar &grounding_line_change, array::Scalar &conservation_error_change, array::Scalar &no_model_mask_change)
Correct the new water thickness based on boundary requirements.
Definition Hydrology.cc:670
const array::Scalar & subglacial_water_thickness() const
Return the effective thickness of the transportable basal water layer.
Definition Hydrology.cc:507
array::Scalar m_basal_melt_rate
Definition Hydrology.hh:182
virtual void init_impl(const array::Scalar &W_till, const array::Scalar &W, const array::Scalar &P)
Definition Hydrology.cc:388
array::Scalar1 m_W
effective thickness of transportable basal water
Definition Hydrology.hh:173
virtual void restart_impl(const File &input_file, int record)
Definition Hydrology.cc:371
array::Scalar m_input_change
Definition Hydrology.hh:193
array::Scalar m_Wtill
effective thickness of basal water stored in till
Definition Hydrology.hh:170
array::Scalar m_Pover
overburden pressure
Definition Hydrology.hh:176
array::Scalar m_grounded_margin_change
Definition Hydrology.hh:191
virtual void write_state_impl(const OutputFile &output) const
The default (empty implementation).
Definition Hydrology.cc:472
array::Scalar m_no_model_mask_change
Definition Hydrology.hh:194
array::Scalar m_conservation_error_change
Definition Hydrology.hh:190
The PISM subglacial hydrology model interface.
Definition Hydrology.hh:109
const Geometry * geometry
Definition Hydrology.hh:37
const array::Scalar1 * no_model_mask
Definition Hydrology.hh:35
virtual void init_impl(const array::Scalar &W_till, const array::Scalar &W, const array::Scalar &P)
Definition Routing.cc:349
array::Staggered1 m_Qstag_average
Definition Routing.hh:115
array::Scalar1 m_R
Definition Routing.hh:130
virtual void write_state_impl(const OutputFile &output) const
The default (empty implementation).
Definition Routing.cc:361
virtual void initialization_message() const
Definition Routing.cc:320
array::Staggered m_Vstag
Definition Routing.hh:118
const array::Scalar & subglacial_water_pressure() const
Definition Routing.cc:368
array::Staggered1 m_Wstag
Definition Routing.hh:121
double max_timestep_W_diff(double KW_max) const
Definition Routing.cc:693
virtual void bootstrap_impl(const File &input_file, const array::Scalar &ice_thickness)
Definition Routing.cc:339
double max_timestep_W_cfl() const
Definition Routing.cc:705
void compute_conductivity(const array::Staggered &W, const array::Scalar &P, const array::Scalar &bed, array::Staggered &result, double &maxKW) const
Compute the nonlinear conductivity at the center of cell edges.
Definition Routing.cc:439
virtual std::map< std::string, Diagnostic::Ptr > spatial_diagnostics_impl() const
Definition Routing.cc:978
array::Staggered1 m_Kstag
Definition Routing.hh:124
void advective_fluxes(const array::Staggered &V, const array::Scalar &W, array::Staggered &result) const
Compute Q = V W at edge-centers (staggered grid) by first-order upwinding.
Definition Routing.cc:673
void W_change_due_to_flow(double dt, const array::Scalar1 &W, const array::Staggered1 &Wstag, const array::Staggered1 &K, const array::Staggered1 &Q, array::Scalar &result)
Definition Routing.cc:765
void update_W(double dt, const array::Scalar &surface_input_rate, const array::Scalar &basal_melt_rate, const array::Scalar1 &W, const array::Staggered1 &Wstag, const array::Scalar &Wtill, const array::Scalar &Wtill_new, const array::Staggered1 &K, const array::Staggered1 &Q, array::Scalar &W_new)
The computation of Wnew, called by update().
Definition Routing.cc:802
virtual void update_impl(double t, double dt, const Inputs &inputs)
Update the model state variables W and Wtill by applying the subglacial hydrology model equations.
Definition Routing.cc:841
void water_thickness_staggered(const array::Scalar &W, const array::CellType1 &mask, array::Staggered &result)
Average the regular grid water thickness to values at the center of cell edges.
Definition Routing.cc:380
void compute_velocity(const array::Staggered &W, const array::Scalar &P, const array::Scalar &bed, const array::Staggered &K, const array::Scalar1 *no_model_mask, array::Staggered &result) const
Get the advection velocity V at the center of cell edges.
Definition Routing.cc:614
const array::Staggered & velocity_staggered() const
Definition Routing.cc:372
array::Staggered1 m_Qstag
Definition Routing.hh:113
virtual std::map< std::string, TSDiagnostic::Ptr > scalar_diagnostics_impl() const
Definition Routing.cc:992
void update_Wtill(double dt, const array::Scalar &Wtill, const array::Scalar &surface_input_rate, const array::Scalar &basal_melt_rate, array::Scalar &Wtill_new)
The computation of Wtillnew, called by update().
Definition Routing.cc:736
array::Scalar m_Wtillnew
Definition Routing.hh:127
Routing(std::shared_ptr< const Grid > g)
Definition Routing.cc:250
array::Scalar m_Wnew
Definition Routing.hh:127
virtual void restart_impl(const File &input_file, int record)
Definition Routing.cc:331
array::Scalar1 m_bottom_surface
Definition Routing.hh:135
virtual std::set< VariableMetadata > state_impl() const
Definition Routing.cc:357
A subglacial hydrology model which assumes water pressure equals overburden pressure.
Definition Routing.hh:81
virtual std::shared_ptr< array::Array > compute_impl() const
Definition Routing.cc:52
Reports the pressure of the transportable water in the subglacial layer.
Definition Routing.cc:43
virtual std::shared_ptr< array::Array > compute_impl() const
Definition Routing.cc:177
Diagnostically reports the staggered-grid components of the velocity of the water in the subglacial l...
Definition Routing.cc:168
virtual std::shared_ptr< array::Array > compute_impl() const
Definition Routing.cc:118
Reports the effective pressure of the transportable water in the subglacial layer,...
Definition Routing.cc:106
std::shared_ptr< array::Array > compute_impl() const
Definition Routing.cc:229
Report hydraulic potential in the subglacial hydrology system.
Definition Routing.cc:221
virtual std::shared_ptr< array::Array > compute_impl() const
Definition Routing.cc:78
Reports the pressure of the transportable water in the subglacial layer as a fraction of the overburd...
Definition Routing.cc:65
virtual std::shared_ptr< array::Array > compute_impl() const
Definition Routing.cc:155
Report the wall melt rate from dissipation of the potential energy of the transportable water.
Definition Routing.cc:142
#define PISM_ERROR_LOCATION
static const double L
Definition exactTestL.cc:40
std::set< VariableMetadata > metadata(std::initializer_list< const Array * > vecs)
Definition Array.cc:1244
void hydraulic_potential(const array::Scalar &W, const array::Scalar &P, const array::Scalar &sea_level, const array::Scalar &bed, const array::Scalar &ice_thickness, array::Scalar &result)
Compute the hydraulic potential.
Definition Routing.cc:190
static double K(double psi_x, double psi_y, double speed, double epsilon)
void wall_melt(const Routing &model, const array::Scalar &bed_elevation, array::Scalar &result)
Compute the wall melt rate which comes from (turbulent) dissipation of flow energy.
Definition Routing.cc:530
void check_bounds(const array::Scalar &W, double W_max)
Definition Hydrology.cc:554
double convert(System::Ptr system, double input, const std::string &spec1, const std::string &spec2)
Convert a quantity from unit1 to unit2.
Definition Units.cc:70
static const double rhow
Definition exactTestP.cc:38
static const double g
Definition exactTestP.cc:36
void GlobalMax(MPI_Comm comm, double *local, double *result, int count)
T clip(T x, T a, T b)
std::map< std::string, Diagnostic::Ptr > DiagnosticList
static const double k
Definition exactTestP.cc:42
void ice_bottom_surface(const Geometry &geometry, array::Scalar &result)
Definition Geometry.cc:241
T combine(const T &a, const T &b)