PISM, A Parallel Ice Sheet Model 2.3.0-79cae578d committed by Constantine Khrulev on 2026-03-22
Loading...
Searching...
No Matches
StressBalance.cc
Go to the documentation of this file.
1// Copyright (C) 2010, 2011, 2012, 2013, 2014, 2015, 2016, 2017, 2018, 2019, 2020, 2021, 2022, 2023, 2024, 2025, 2026 Constantine Khroulev and Ed Bueler
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#include <memory>
19
20#include "pism/stressbalance/StressBalance.hh"
21#include "pism/geometry/Geometry.hh"
22#include "pism/rheology/FlowLaw.hh"
23#include "pism/stressbalance/SSB_Modifier.hh"
24#include "pism/stressbalance/ShallowStressBalance.hh"
25#include "pism/util/Config.hh"
26#include "pism/util/Context.hh"
27#include "pism/util/EnthalpyConverter.hh"
28#include "pism/util/Grid.hh"
29#include "pism/util/Mask.hh"
30#include "pism/util/Profiling.hh"
31#include "pism/util/Time.hh"
32#include "pism/util/array/CellType.hh"
33#include "pism/util/error_handling.hh"
34#include "pism/util/io/SynchronousOutputWriter.hh"
35#include "pism/util/io/io_helpers.hh"
36
37namespace pism {
38namespace stressbalance {
39
41 geometry = NULL;
42 new_bed_elevation = true;
43
44 basal_melt_rate = NULL;
46 fracture_density = NULL;
47 basal_yield_stress = NULL;
48
49 enthalpy = NULL;
50 age = NULL;
51
52 bc_mask = NULL;
53 bc_values = NULL;
54
55 no_model_mask = NULL;
58}
59
60/*!
61 * Save stress balance inputs to a file (for debugging).
62 */
63void Inputs::dump(const char *filename) const {
64 if (geometry == nullptr) {
65 return;
66 }
67
68 auto grid = geometry->ice_thickness.grid();
69 auto ctx = grid->ctx();
70 auto config = ctx->config();
71
72 VariableMetadata mapping{ "mapping", ctx->unit_system() };
73 auto writer = std::make_shared<SynchronousOutputWriter>(ctx->com(), *config);
74 writer->initialize({}, true);
75
76 OutputFile output(writer, filename);
77
78 io::write_config(*config, "pism_config", output);
79
80 auto time = ctx->time();
81 output.define_variable(time->metadata());
82 output.append_time(time->current());
83
84 const array::Array *geom[] = { &geometry->latitude,
93
94 // define
95 for (const auto * vec : geom) {
96 for (const auto &var : vec->all_metadata()) {
97 output.define_variable(var);
98 }
99 }
100 // write
101 for (const auto * vec : geom) {
102 vec->write(output);
103 }
104
105 const array::Array *optional[] = { basal_melt_rate,
109 enthalpy,
110 age,
111 bc_mask,
112 bc_values,
116 // define
117 for (const auto * vec : optional) {
118 for (const auto &var : vec->all_metadata()) {
119 output.define_variable(var);
120 }
121 }
122 // write
123 for (const auto * vec : optional) {
124 vec->write(output);
125 }
126}
127
128StressBalance::StressBalance(std::shared_ptr<const Grid> g,
129 std::shared_ptr<ShallowStressBalance> sb,
130 std::shared_ptr<SSB_Modifier> ssb_mod)
131 : Component(g),
132 m_w(m_grid, "wvel_rel", array::WITHOUT_GHOSTS, m_grid->z()),
133 m_strain_heating(m_grid, "strain_heating", array::WITHOUT_GHOSTS, m_grid->z()),
134 m_shallow_stress_balance(sb),
135 m_modifier(ssb_mod) {
136
137 m_w.metadata(0)
138 .long_name("vertical velocity of ice, relative to base of ice directly below")
139 .units("m s^-1")
140 .output_units("m year^-1")
141 .set_time_dependent(true);
142
144 .long_name("rate of strain heating in ice (dissipation heating)")
145 .units("W m^-3");
146}
147
150
151//! \brief Initialize the StressBalance object.
154 m_modifier->init();
155}
156
157//! \brief Performs the shallow stress balance computation.
158void StressBalance::update(const Inputs &inputs, bool full_update) {
159
160 try {
161 profiling().begin("stress_balance.shallow");
162 m_shallow_stress_balance->update(inputs, full_update);
163 profiling().end("stress_balance.shallow");
164
165 profiling().begin("stress_balance.modifier");
166 m_modifier->update(m_shallow_stress_balance->velocity(),
167 inputs, full_update);
168 profiling().end("stress_balance.modifier");
169
170 if (full_update) {
171 const array::Array3D &u = m_modifier->velocity_u();
172 const array::Array3D &v = m_modifier->velocity_v();
173
174 profiling().begin("stress_balance.strain_heat");
176 profiling().end("stress_balance.strain_heat");
177
178 profiling().begin("stress_balance.vertical_velocity");
180 u, v, inputs.basal_melt_rate, m_w);
181 profiling().end("stress_balance.vertical_velocity");
182
184 inputs.geometry->cell_type,
185 inputs.no_model_mask,
186 u, v, m_w);
187 }
188
190 inputs.geometry->cell_type,
191 inputs.no_model_mask,
192 m_shallow_stress_balance->velocity());
193 }
194 catch (RuntimeError &e) {
195 e.add_context("updating the stress balance");
196 throw;
197 }
198}
199
203
207
211
213 return m_modifier->diffusive_flux();
214}
215
217 return m_modifier->max_diffusivity();
218}
219
221 return m_modifier->velocity_u();
222}
223
225 return m_modifier->velocity_v();
226}
227
229 return m_w;
230}
231
233 return m_shallow_stress_balance->basal_frictional_heating();
234}
235
239
240//! Compute vertical velocity using incompressibility of the ice.
241/*!
242The vertical velocity \f$w(x,y,z,t)\f$ is the velocity *relative to the
243location of the base of the ice column*. That is, the vertical velocity
244computed here is identified as \f$\tilde w(x,y,s,t)\f$ in the page
245[]@ref vertchange.
246
247Thus \f$w<0\f$ here means that that
248that part of the ice is getting closer to the base of the ice, and so on.
249The slope of the bed (i.e. relative to the geoid) and/or the motion of the
250bed (i.e. from bed deformation) do not affect the vertical velocity.
251
252In fact the following statement is exactly true if the basal melt rate is zero:
253the vertical velocity at a point in the ice is positive (negative) if and only
254if the average horizontal divergence of the horizontal velocity, in the portion
255of the ice column below that point, is negative (positive).
256In particular, because \f$z=0\f$ is the location of the base of the ice
257always, the only way to have \f$w(x,y,0,t) \ne 0\f$ is to have a basal melt
258rate.
259
260Incompressibility itself says
261 \f[ \nabla\cdot\mathbf{U} + \frac{\partial w}{\partial z} = 0. \f]
262This is immediately equivalent to the integral
263 \f[ w(x,y,z,t) = - \int_{b(x,y,t)}^{z} \nabla\cdot\mathbf{U}\,d\zeta
264 + w_b(x,y,t). \f]
265Here the value \f$w_b(x,y,t)\f$ is either zero or the negative of the basal melt rate
266according to the value of the flag `geometry.update.use_basal_melt_rate`.
267
268The vertical integral is computed by the trapezoid rule.
269 */
271 const array::Array3D &u,
272 const array::Array3D &v,
273 const array::Scalar *basal_melt_rate,
274 array::Array3D &result) {
275
276 const bool use_upstream_fd = m_config->get_string("stress_balance.vertical_velocity_approximation") == "upstream";
277
278 array::AccessScope list{&u, &v, &mask, &result};
279
280 if (basal_melt_rate) {
281 list.add(*basal_melt_rate);
282 }
283
284 const std::vector<double> &z = m_grid->z();
285 const unsigned int Mz = m_grid->Mz();
286
287 const double
288 dx = m_grid->dx(),
289 dy = m_grid->dy();
290
291 std::vector<double> u_x_plus_v_y(Mz);
292
293 for (auto p : m_grid->points()) {
294 const int i = p.i(), j = p.j();
295
296 double *w_ij = result.get_column(i,j);
297
298 const double
299 *u_w = u.get_column(i-1,j),
300 *u_ij = u.get_column(i,j),
301 *u_e = u.get_column(i+1,j);
302 const double
303 *v_s = v.get_column(i,j-1),
304 *v_ij = v.get_column(i,j),
305 *v_n = v.get_column(i,j+1);
306
307 double
308 west = 1.0,
309 east = 1.0,
310 south = 1.0,
311 north = 1.0;
312 double
313 D_x = 0, // 1/(dx), 1/(2dx), or 0
314 D_y = 0; // 1/(dy), 1/(2dy), or 0
315
316 // Switch between second-order centered differences in the interior and
317 // first-order one-sided differences at ice margins.
318
319 // x-derivative
320 {
321 // use basal velocity to determine FD direction ("upwind" when it's clear, centered when it's
322 // not)
323 if (use_upstream_fd) {
324 const double
325 uw = 0.5 * (u_w[0] + u_ij[0]),
326 ue = 0.5 * (u_ij[0] + u_e[0]);
327
328 if (uw > 0.0 and ue >= 0.0) {
329 west = 1.0;
330 east = 0.0;
331 } else if (uw <= 0.0 and ue < 0.0) {
332 west = 0.0;
333 east = 1.0;
334 } else {
335 west = 1.0;
336 east = 1.0;
337 }
338 }
339
340 if ((mask.icy(i,j) and mask.ice_free(i+1,j)) or (mask.ice_free(i,j) and mask.icy(i+1,j))) {
341 east = 0;
342 }
343 if ((mask.icy(i,j) and mask.ice_free(i-1,j)) or (mask.ice_free(i,j) and mask.icy(i-1,j))) {
344 west = 0;
345 }
346
347 if (east + west > 0) {
348 D_x = 1.0 / (dx * (east + west));
349 } else {
350 D_x = 0.0;
351 }
352 }
353
354 // y-derivative
355 {
356 // use basal velocity to determine FD direction ("upwind" when it's clear, centered when it's
357 // not)
358 if (use_upstream_fd) {
359 const double
360 vs = 0.5 * (v_s[0] + v_ij[0]),
361 vn = 0.5 * (v_ij[0] + v_n[0]);
362
363 if (vs > 0.0 and vn >= 0.0) {
364 south = 1.0;
365 north = 0.0;
366 } else if (vs <= 0.0 and vn < 0.0) {
367 south = 0.0;
368 north = 1.0;
369 } else {
370 south = 1.0;
371 north = 1.0;
372 }
373 }
374
375 if ((mask.icy(i,j) and mask.ice_free(i,j+1)) or (mask.ice_free(i,j) and mask.icy(i,j+1))) {
376 north = 0;
377 }
378 if ((mask.icy(i,j) and mask.ice_free(i,j-1)) or (mask.ice_free(i,j) and mask.icy(i,j-1))) {
379 south = 0;
380 }
381
382 if (north + south > 0) {
383 D_y = 1.0 / (dy * (north + south));
384 } else {
385 D_y = 0.0;
386 }
387 }
388
389 // compute u_x + v_y using a vectorizable loop
390 for (unsigned int k = 0; k < Mz; ++k) {
391 double
392 u_x = D_x * (west * (u_ij[k] - u_w[k]) + east * (u_e[k] - u_ij[k])),
393 v_y = D_y * (south * (v_ij[k] - v_s[k]) + north * (v_n[k] - v_ij[k]));
394 u_x_plus_v_y[k] = u_x + v_y;
395 }
396
397 // at the base: include the basal melt rate
398 if (basal_melt_rate != NULL) {
399 w_ij[0] = - (*basal_melt_rate)(i,j);
400 } else {
401 w_ij[0] = 0.0;
402 }
403
404 // within the ice and above:
405 for (unsigned int k = 1; k < Mz; ++k) {
406 const double dz = z[k] - z[k-1];
407
408 w_ij[k] = w_ij[k - 1] - (0.5 * dz) * (u_x_plus_v_y[k] + u_x_plus_v_y[k - 1]);
409 }
410 }
411}
412
413/**
414 * This function computes \f$D^2\f$ defined by
415 *
416 * \f[ 2D^2 = D_{ij} D_{ij}\f]
417 * or
418 * \f[
419 * D^2 = \frac{1}{2}\,\left(\frac{1}{2}\,(v_{z})^2 + (v_{y} + u_{x})^2 +
420 * (v_{y})^2 + \frac{1}{2}\,(v_{x} + u_{y})^2 + \frac{1}{2}\,(u_{z})^2 +
421 * (u_{x})^2\right)
422 * \f]
423 *
424 * (note the use of the summation convention). Here \f$D_{ij}\f$ is the
425 * strain rate tensor. See
426 * StressBalance::compute_volumetric_strain_heating() for details.
427 *
428 * @param u_x,u_y,u_z partial derivatives of \f$u\f$, the x-component of the ice velocity
429 * @param v_x,v_y,v_z partial derivatives of \f$v\f$, the y-component of the ice velocity
430 *
431 * @return \f$D^2\f$, where \f$D\f$ is defined above.
432 */
433static inline double D2(double u_x, double u_y, double u_z, double v_x, double v_y, double v_z) {
434 return 0.5 * (PetscSqr(u_x + v_y) + u_x*u_x + v_y*v_y + 0.5 * (PetscSqr(u_y + v_x) + u_z*u_z + v_z*v_z));
435}
436
437/**
438 \brief Computes the volumetric strain heating using horizontal
439 velocity.
440
441 Following the notation used in [\ref BBssasliding], let \f$u\f$ be a
442 three-dimensional *vector* velocity field. Then the strain rate
443 tensor \f$D_{ij}\f$ is defined by
444
445 \f[ D_{ij} = \frac 12 \left(\diff{u_{i}}{x_{j}} + \diff{u_{j}}{x_{i}} \right), \f]
446
447 Where \f$i\f$ and \f$j\f$ range from \f$1\f$ to \f$3\f$.
448
449 The flow law in the viscosity form states
450
451 \f[ \tau_{ij} = 2 \eta D_{ij}, \f]
452
453 and the nonlinear ice viscosity satisfies
454
455 \f[ 2 \eta = B(T) D^{(1/n) - 1}. \f]
456
457 Here \f$D^{2}\f$ is defined by \f$2D^{2} = D_{ij}D_{ij}\f$ (using the
458 summation convention) and \f$B(T) = A(T)^{-1/n}\f$ is the ice hardness.
459
460 Now the volumetric strain heating is
461
462 \f[ \Sigma = \sum_{i,j=1}^{3}D_{ij}\tau_{ij} = 2 B(T) D^{(1/n) + 1}. \f]
463
464 We use an *approximation* of \f$D_{ij}\f$ common in shallow ice models:
465
466 - we assume that horizontal derivatives of the vertical velocity are
467 much smaller than \f$z\f$ derivatives horizontal velocity
468 components \f$u\f$ and \f$v\f$. (We drop \f$w_x\f$ and \f$w_y\f$
469 terms in \f$D_{ij}\f$.)
470
471 - we use the incompressibility of ice to approximate \f$w_z\f$:
472
473 \f[ w_z = - (u_x + v_y). \f]
474
475 Requires ghosts of `u` and `v` velocity components and uses the fact
476 that `u` and `v` above the ice are filled using constant
477 extrapolation.
478
479 Resulting field does not have ghosts.
480
481 Below is the *Maxima* code that produces the expression evaluated by D2().
482
483 derivabbrev : true;
484 U : [u, v, w]; X : [x, y, z]; depends(U, X);
485 gradef(w, x, 0); gradef(w, y, 0);
486 gradef(w, z, -(diff(u, x) + diff(v, y)));
487 d[i,j] := 1/2 * (diff(U[i], X[j]) + diff(U[j], X[i]));
488 D : genmatrix(d, 3, 3), ratsimp, factor;
489 tex('D = D);
490 tex('D^2 = 1/2 * mat_trace(D . D));
491
492 @return 0 on success
493 */
495 PetscErrorCode ierr;
496
497 const auto &flow_law = *m_shallow_stress_balance->flow_law();
498 auto EC = m_shallow_stress_balance->enthalpy_converter();
499
500 const array::Array3D
501 &u = m_modifier->velocity_u(),
502 &v = m_modifier->velocity_v();
503
504 const array::Scalar &thickness = inputs.geometry->ice_thickness;
505 const array::Array3D *enthalpy = inputs.enthalpy;
506
507 const auto &mask = inputs.geometry->cell_type;
508
509 double
510 enhancement_factor = m_shallow_stress_balance->flow_enhancement_factor(),
511 n = flow_law.exponent(),
512 exponent = 0.5 * (1.0 / n + 1.0),
513 e_to_a_power = pow(enhancement_factor,-1.0/n);
514
515 array::AccessScope list{&mask, enthalpy, &m_strain_heating, &thickness, &u, &v};
516
517 const std::vector<double> &z = m_grid->z();
518 const unsigned int Mz = m_grid->Mz();
519 std::vector<double> depth(Mz), pressure(Mz), hardness(Mz);
520
521 ParallelSection loop(m_grid->com);
522 try {
523 for (auto p : m_grid->points()) {
524 const int i = p.i(), j = p.j();
525
526 double H = thickness(i, j);
527 int ks = m_grid->kBelowHeight(H);
528 const double
529 *u_ij, *u_w, *u_n, *u_e, *u_s,
530 *v_ij, *v_w, *v_n, *v_e, *v_s;
531 double *Sigma;
532 const double *E_ij;
533
534 double west = 1, east = 1, south = 1, north = 1,
535 D_x = 0, // 1/(dx), 1/(2dx), or 0
536 D_y = 0; // 1/(dy), 1/(2dy), or 0
537
538 // x-derivative
539 {
540 if ((mask.icy(i,j) and mask.ice_free(i+1,j)) or (mask.ice_free(i,j) and mask.icy(i+1,j))) {
541 east = 0;
542 }
543 if ((mask.icy(i,j) and mask.ice_free(i-1,j)) or (mask.ice_free(i,j) and mask.icy(i-1,j))) {
544 west = 0;
545 }
546
547 if (east + west > 0) {
548 D_x = 1.0 / (m_grid->dx() * (east + west));
549 } else {
550 D_x = 0.0;
551 }
552 }
553
554 // y-derivative
555 {
556 if ((mask.icy(i,j) and mask.ice_free(i,j+1)) or (mask.ice_free(i,j) and mask.icy(i,j+1))) {
557 north = 0;
558 }
559 if ((mask.icy(i,j) and mask.ice_free(i,j-1)) or (mask.ice_free(i,j) and mask.icy(i,j-1))) {
560 south = 0;
561 }
562
563 if (north + south > 0) {
564 D_y = 1.0 / (m_grid->dy() * (north + south));
565 } else {
566 D_y = 0.0;
567 }
568 }
569
570 u_ij = u.get_column(i, j);
571 u_w = u.get_column(i - 1, j);
572 u_e = u.get_column(i + 1, j);
573 u_s = u.get_column(i, j - 1);
574 u_n = u.get_column(i, j + 1);
575
576 v_ij = v.get_column(i, j);
577 v_w = v.get_column(i - 1, j);
578 v_e = v.get_column(i + 1, j);
579 v_s = v.get_column(i, j - 1);
580 v_n = v.get_column(i, j + 1);
581
582 E_ij = enthalpy->get_column(i, j);
583 Sigma = m_strain_heating.get_column(i, j);
584
585 for (int k = 0; k <= ks; ++k) {
586 depth[k] = H - z[k];
587 }
588
589 // pressure added by the ice (i.e. pressure difference between the
590 // current level and the top of the column)
591 EC->pressure(depth, ks, pressure); // FIXME issue #15
592
593 flow_law.hardness_n(E_ij, pressure.data(), ks + 1, hardness.data());
594
595 for (int k = 0; k <= ks; ++k) {
596 double dz;
597
598 double u_z = 0.0, v_z = 0.0,
599 u_x = D_x * (west * (u_ij[k] - u_w[k]) + east * (u_e[k] - u_ij[k])),
600 u_y = D_y * (south * (u_ij[k] - u_s[k]) + north * (u_n[k] - u_ij[k])),
601 v_x = D_x * (west * (v_ij[k] - v_w[k]) + east * (v_e[k] - v_ij[k])),
602 v_y = D_y * (south * (v_ij[k] - v_s[k]) + north * (v_n[k] - v_ij[k]));
603
604 if (k > 0) {
605 dz = z[k+1] - z[k-1];
606 u_z = (u_ij[k+1] - u_ij[k-1]) / dz;
607 v_z = (v_ij[k+1] - v_ij[k-1]) / dz;
608 } else {
609 // use one-sided differences for u_z and v_z on the bottom level
610 dz = z[1] - z[0];
611 u_z = (u_ij[1] - u_ij[0]) / dz;
612 v_z = (v_ij[1] - v_ij[0]) / dz;
613 }
614
615 Sigma[k] = 2.0 * e_to_a_power * hardness[k] * pow(D2(u_x, u_y, u_z, v_x, v_y, v_z), exponent);
616 } // k-loop
617
618 int remaining_levels = Mz - (ks + 1);
619 if (remaining_levels > 0) {
620#if PETSC_VERSION_LT(3, 12, 0)
621 ierr = PetscMemzero(&Sigma[ks+1],
622 remaining_levels*sizeof(double));
623#else
624 ierr = PetscArrayzero(&Sigma[ks+1], remaining_levels);
625#endif
626 PISM_CHK(ierr, "PetscMemzero");
627 }
628 }
629 } catch (...) {
630 loop.failed();
631 }
632 loop.check();
633}
634
635std::string StressBalance::stdout_report() const {
636 return m_shallow_stress_balance->stdout_report() + m_modifier->stdout_report();
637}
638
642
644 return m_modifier.get();
645}
646
647std::set<VariableMetadata> StressBalance::state_impl() const {
649 auto modifier = m_modifier->state();
650
652}
653
655 m_shallow_stress_balance->write_state(output);
656 m_modifier->write_state(output);
657}
658
659//! \brief Compute eigenvalues of the horizontal, vertically-integrated strain rate tensor.
660/*!
661Calculates all components \f$D_{xx}, D_{yy}, D_{xy}=D_{yx}\f$ of the
662vertically-averaged strain rate tensor \f$D\f$ [\ref SchoofStream]. Then computes
663the eigenvalues `result(i,j,0)` = (maximum eigenvalue), `result(i,j,1)` = (minimum
664eigenvalue). Uses the provided thickness to make decisions (PIK) about computing
665strain rates near calving front.
666
667Note that `result(i,j,0)` >= `result(i,j,1)`, but there is no necessary relation between
668the magnitudes, and either principal strain rate could be negative or positive.
669
670Result can be used in a calving law, for example in eigencalving (PIK).
671
672Note: strain rates will be derived from SSA velocities, using ghosts when
673necessary. Both implementations (SSAFD and SSAFEM) call
674update_ghosts() to ensure that ghost values are up to date.
675 */
677 const array::CellType1 &mask,
679
680 using mask::ice_free;
681
682 auto grid = result.grid();
683 double dx = grid->dx();
684 double dy = grid->dy();
685
686 array::AccessScope list{&V, &mask, &result};
687
688 for (auto p : grid->points()) {
689 const int i = p.i(), j = p.j();
690
691 if (mask.ice_free(i,j)) {
692 result(i, j).eigen1 = 0.0;
693 result(i, j).eigen2 = 0.0;
694 continue;
695 }
696
697 auto m = mask.star(i,j);
698 auto U = V.star(i,j);
699
700 // strain in units s^-1
701 double u_x = 0, u_y = 0, v_x = 0, v_y = 0,
702 east = 1, west = 1, south = 1, north = 1;
703
704 // Computes u_x using second-order centered finite differences written as
705 // weighted sums of first-order one-sided finite differences.
706 //
707 // Given the cell layout
708 // *----n----*
709 // | |
710 // | |
711 // w e
712 // | |
713 // | |
714 // *----s----*
715 // east == 0 if the east neighbor of the current cell is ice-free. In
716 // this case we use the left- (west-) sided difference.
717 //
718 // If both neighbors in the east-west (x) direction are ice-free the
719 // x-derivative is set to zero (see u_x, v_x initialization above).
720 //
721 // Similarly in other directions.
722 if (ice_free(m.e)) {
723 east = 0;
724 }
725 if (ice_free(m.w)) {
726 west = 0;
727 }
728 if (ice_free(m.n)) {
729 north = 0;
730 }
731 if (ice_free(m.s)) {
732 south = 0;
733 }
734
735 if (west + east > 0) {
736 u_x = 1.0 / (dx * (west + east)) * (west * (U.c.u - U[West].u) + east * (U[East].u - U.c.u));
737 v_x = 1.0 / (dx * (west + east)) * (west * (U.c.v - U[West].v) + east * (U[East].v - U.c.v));
738 }
739
740 if (south + north > 0) {
741 u_y = 1.0 / (dy * (south + north)) * (south * (U.c.u - U[South].u) + north * (U[North].u - U.c.u));
742 v_y = 1.0 / (dy * (south + north)) * (south * (U.c.v - U[South].v) + north * (U[North].v - U.c.v));
743 }
744
745 const double A = 0.5 * (u_x + v_y), // A = (1/2) trace(D)
746 B = 0.5 * (u_x - v_y),
747 Dxy = 0.5 * (v_x + u_y), // B^2 = A^2 - u_x v_y
748 q = sqrt(B*B + Dxy*Dxy);
749 result(i, j).eigen1 = A + q;
750 result(i, j).eigen2 = A - q; // q >= 0 so e1 >= e2
751
752 }
753}
754
755//! @brief Compute 2D deviatoric stresses.
757 const array::Vector1 &velocity,
758 const array::Scalar &hardness,
759 const array::CellType1 &cell_type,
761
762 using mask::ice_free;
763
764 auto grid = result.grid();
765
766 const double
767 dx = grid->dx(),
768 dy = grid->dy();
769
770 array::AccessScope list{&velocity, &hardness, &result, &cell_type};
771
772 for (auto p : grid->points()) {
773 const int i = p.i(), j = p.j();
774
775 if (cell_type.ice_free(i, j)) {
776 result(i,j).xx = 0.0;
777 result(i,j).yy = 0.0;
778 result(i,j).xy = 0.0;
779 continue;
780 }
781
782 auto m = cell_type.star(i,j);
783 auto U = velocity.star(i,j);
784
785 // strain in units s^-1
786 double u_x = 0, u_y = 0, v_x = 0, v_y = 0,
787 east = 1, west = 1, south = 1, north = 1;
788
789 // Computes u_x using second-order centered finite differences written as
790 // weighted sums of first-order one-sided finite differences.
791 //
792 // Given the cell layout
793 // *----n----*
794 // | |
795 // | |
796 // w e
797 // | |
798 // | |
799 // *----s----*
800 // east == 0 if the east neighbor of the current cell is ice-free. In
801 // this case we use the left- (west-) sided difference.
802 //
803 // If both neighbors in the east-west (x) direction are ice-free the
804 // x-derivative is set to zero (see u_x, v_x initialization above).
805 //
806 // Similarly in y-direction.
807 if (ice_free(m.e)) {
808 east = 0;
809 }
810 if (ice_free(m.w)) {
811 west = 0;
812 }
813 if (ice_free(m.n)) {
814 north = 0;
815 }
816 if (ice_free(m.s)) {
817 south = 0;
818 }
819
820 if (west + east > 0) {
821 u_x = 1.0 / (dx * (west + east)) * (west * (U.c.u - U[West].u) + east * (U[East].u - U.c.u));
822 v_x = 1.0 / (dx * (west + east)) * (west * (U.c.v - U[West].v) + east * (U[East].v - U.c.v));
823 }
824
825 if (south + north > 0) {
826 u_y = 1.0 / (dy * (south + north)) * (south * (U.c.u - U[South].u) + north * (U[North].u - U.c.u));
827 v_y = 1.0 / (dy * (south + north)) * (south * (U.c.v - U[South].v) + north * (U[North].v - U.c.v));
828 }
829
830 double nu = 0.0;
831 flow_law.effective_viscosity(hardness(i, j),
832 secondInvariant_2D({u_x, v_x}, {u_y, v_y}),
833 &nu, NULL);
834
835 //get deviatoric stresses
836 result(i,j).xx = 2.0*nu*u_x;
837 result(i,j).yy = 2.0*nu*v_y;
838 result(i,j).xy = nu*(u_y+v_x);
839 }
840}
841
842
843} // end of namespace stressbalance
844} // end of namespace pism
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
const Profiling & profiling() const
Definition Component.cc:115
std::set< VariableMetadata > state() const
Definition Component.cc:124
A class defining a common interface for most PISM sub-models.
Definition Component.hh:118
array::Scalar1 sea_level_elevation
Definition Geometry.hh:48
array::Scalar cell_grounded_fraction
Definition Geometry.hh:56
array::Scalar2 ice_surface_elevation
Definition Geometry.hh:57
array::Scalar1 ice_area_specific_volume
Definition Geometry.hh:52
array::CellType2 cell_type
Definition Geometry.hh:55
array::Scalar2 ice_thickness
Definition Geometry.hh:51
array::Scalar longitude
Definition Geometry.hh:42
array::Scalar2 bed_elevation
Definition Geometry.hh:47
array::Scalar latitude
Definition Geometry.hh:41
void append_time(double time_seconds) const
Definition OutputFile.cc:41
void define_variable(const VariableMetadata &variable) const
Definition OutputFile.cc:31
void failed()
Indicates a failure of a parallel section.
void begin(const char *name) const
Definition Profiling.cc:75
void end(const char *name) const
Definition Profiling.cc:91
void add_context(const std::string &message)
Add a message providing some context. This way we can (sort of) get a stack trace even though C++ exc...
VariableMetadata & long_name(const std::string &input)
VariableMetadata & units(const std::string &input)
VariableMetadata & output_units(const std::string &input)
std::shared_ptr< units::System > unit_system() const
VariableMetadata & set_time_dependent(bool flag)
Makes sure that we call begin_access() and end_access() for all accessed array::Arrays.
Definition Array.hh:66
stencils::Star< T > star(int i, int j) const
Definition Array2D.hh:107
A storage vector combining related fields in a struct.
Definition Array2D.hh:32
double * get_column(int i, int j)
Definition Array3D.cc:127
A virtual class collecting methods common to ice and bedrock 3D fields.
Definition Array3D.hh:33
std::shared_ptr< const Grid > grid() const
Definition Array.cc:134
void add(double alpha, const Array &x)
Result: v <- v + alpha * x. Calls VecAXPY.
Definition Array.cc:209
VariableMetadata & metadata(unsigned int N=0)
Returns a reference to the VariableMetadata object containing metadata for the compoment N.
Definition Array.cc:467
Abstract class for reading, writing, allocating, and accessing a DA-based PETSc Vec (2D and 3D fields...
Definition Array.hh:209
bool ice_free(int i, int j) const
Definition CellType.hh:54
bool icy(int i, int j) const
Definition CellType.hh:42
A class for storing and accessing internal staggered-grid 2D fields. Uses dof=2 storage....
Definition Staggered.hh:37
void effective_viscosity(double hardness, double gamma, double *nu, double *dnu) const
Computes the regularized effective viscosity and its derivative with respect to the second invariant ...
Definition FlowLaw.cc:170
const array::Scalar1 * fracture_density
const array::Scalar * basal_yield_stress
const array::Scalar * bc_mask
const array::Array3D * age
void dump(const char *filename) const
const array::Scalar * no_model_ice_thickness
const array::Scalar2 * no_model_surface_elevation
const array::Scalar * water_column_pressure
const array::Array3D * enthalpy
const array::Scalar2 * no_model_mask
const array::Scalar * basal_melt_rate
const array::Vector * bc_values
Shallow stress balance modifier (such as the non-sliding SIA).
Shallow stress balance (such as the SSA).
std::shared_ptr< SSB_Modifier > m_modifier
const array::Array3D & velocity_w() const
virtual void write_state_impl(const OutputFile &output) const
The default (empty implementation).
const SSB_Modifier * modifier() const
Returns a pointer to a stress balance modifier implementation.
const array::Array3D & velocity_u() const
Get components of the the 3D velocity field.
virtual void compute_volumetric_strain_heating(const Inputs &inputs)
Computes the volumetric strain heating using horizontal velocity.
std::shared_ptr< ShallowStressBalance > m_shallow_stress_balance
double max_diffusivity() const
Get the max diffusivity (for the adaptive time-stepping).
virtual void compute_vertical_velocity(const array::CellType1 &mask, const array::Array3D &u, const array::Array3D &v, const array::Scalar *bmr, array::Array3D &result)
Compute vertical velocity using incompressibility of the ice.
void init()
Initialize the StressBalance object.
const array::Vector & advective_velocity() const
Get the thickness-advective (SSA) 2D velocity.
const array::Staggered & diffusive_flux() const
Get the diffusive (SIA) vertically-averaged flux on the staggered grid.
const ShallowStressBalance * shallow() const
Returns a pointer to a shallow stress balance solver implementation.
std::string stdout_report() const
Produce a report string for the standard output.
void update(const Inputs &inputs, bool full_update)
Update all the fields if (full_update), only update diffusive flux and max. diffusivity otherwise.
const array::Scalar & basal_frictional_heating() const
Get the basal frictional heating.
const array::Array3D & velocity_v() const
const array::Array3D & volumetric_strain_heating() const
StressBalance(std::shared_ptr< const Grid > g, std::shared_ptr< ShallowStressBalance > sb, std::shared_ptr< SSB_Modifier > ssb_mod)
virtual std::set< VariableMetadata > state_impl() const
#define PISM_CHK(errcode, name)
#define n
Definition exactTestM.c:37
void write_config(const Config &config, const std::string &variable_name, const OutputFile &file)
bool ice_free(int M)
Ice-free cell (grounded or ocean).
Definition Mask.hh:58
void compute_2D_principal_strain_rates(const array::Vector1 &V, const array::CellType1 &mask, array::Array2D< PrincipalStrainRates > &result)
Compute eigenvalues of the horizontal, vertically-integrated strain rate tensor.
static double D2(double u_x, double u_y, double u_z, double v_x, double v_y, double v_z)
void compute_2D_stresses(const rheology::FlowLaw &flow_law, const array::Vector1 &velocity, const array::Scalar &hardness, const array::CellType1 &cell_type, array::Array2D< DeviatoricStresses > &result)
Compute 2D deviatoric stresses.
static double secondInvariant_2D(const Vector2d &U_x, const Vector2d &U_y)
Definition FlowLaw.hh:45
static const double g
Definition exactTestP.cc:36
CFLData max_timestep_cfl_2d(const array::Scalar &ice_thickness, const array::CellType &cell_type, const array::Scalar1 *no_model_mask, const array::Vector &velocity)
Compute the CFL constant associated to first-order upwinding for the sliding contribution to mass con...
static const double k
Definition exactTestP.cc:42
@ North
Definition stencils.hh:24
@ East
Definition stencils.hh:24
@ South
Definition stencils.hh:24
@ West
Definition stencils.hh:24
T combine(const T &a, const T &b)
CFLData max_timestep_cfl_3d(const array::Scalar &ice_thickness, const array::CellType &cell_type, const array::Scalar1 *no_model_mask, const array::Array3D &u3, const array::Array3D &v3, const array::Array3D &w3)
Compute the maximum velocities for time-stepping and reporting to user.