PISM, A Parallel Ice Sheet Model 2.3.0-79cae578d committed by Constantine Khrulev on 2026-03-22
Loading...
Searching...
No Matches
SSAFEM.cc
Go to the documentation of this file.
1// Copyright (C) 2009--2025 Jed Brown and Ed Bueler and Constantine Khroulev and David Maxwell
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 "pism/util/Grid.hh"
20#include "pism/stressbalance/ssa/SSAFEM.hh"
21#include "pism/util/fem/Quadrature.hh"
22#include "pism/util/fem/DirichletData.hh"
23#include "pism/util/Mask.hh"
24#include "pism/basalstrength/basal_resistance.hh"
25#include "pism/rheology/FlowLaw.hh"
26#include "pism/util/pism_options.hh"
27#include "pism/util/error_handling.hh"
28#include "pism/util/Vars.hh"
29#include "pism/stressbalance/StressBalance.hh"
30#include "pism/geometry/Geometry.hh"
31
32#include "pism/util/node_types.hh"
33#include "pism/util/pism_utilities.hh" // average_water_column_pressure()
34#include "pism/util/Interpolation1D.hh"
35#include "pism/util/petscwrappers/DM.hh"
36#include "pism/util/petscwrappers/Vec.hh"
37#include "pism/util/petscwrappers/Viewer.hh"
38#include "pism/util/Logger.hh"
39
40namespace pism {
41namespace stressbalance {
42
43/** The Q1 finite element SSA solver.
44 *
45 *
46 *
47 */
48SSAFEM::SSAFEM(std::shared_ptr<const Grid> grid)
49 : SSA(grid),
50 m_bc_mask(grid, "bc_mask"),
51 m_bc_values(grid, "_bc"),
52 m_gc(*m_config),
53 m_coefficients(grid, "ssa_coefficients", array::WITH_GHOSTS, 1),
54 m_node_type(m_grid, "node_type"),
55 m_boundary_integral(m_grid, "boundary_integral"),
56 m_element_index(*grid),
57 m_q1_element(*grid, fem::Q1Quadrature4()) {
60
61 const double ice_density = m_config->get_number("constants.ice.density");
62 m_alpha = 1 - ice_density / m_config->get_number("constants.sea_water.density");
63 m_rho_g = ice_density * m_config->get_number("constants.standard_gravity");
64
65 m_driving_stress_x = NULL;
66 m_driving_stress_y = NULL;
67
68 PetscErrorCode ierr;
69
70 m_dirichletScale = 1.0;
71 m_beta_ice_free_bedrock = m_config->get_number("basal_resistance.beta_ice_free_bedrock");
72
73 ierr = SNESCreate(m_grid->com, m_snes.rawptr());
74 PISM_CHK(ierr, "SNESCreate");
75
76 // Set the SNES callbacks to call into our compute_local_function and compute_local_jacobian.
78 m_callback_data.ssa = this;
79
80 auto dm = m_velocity_global.dm();
81
82 ierr = DMDASNESSetFunctionLocal(*dm, INSERT_VALUES,
83#if PETSC_VERSION_LT(3,21,0)
84 (DMDASNESFunction)function_callback,
85#else
86 (DMDASNESFunctionFn*)function_callback,
87#endif
89 PISM_CHK(ierr, "DMDASNESSetFunctionLocal");
90
91 ierr = DMDASNESSetJacobianLocal(*dm,
92#if PETSC_VERSION_LT(3,21,0)
93 (DMDASNESJacobian)jacobian_callback,
94#else
95 (DMDASNESJacobianFn*)jacobian_callback,
96#endif
98 PISM_CHK(ierr, "DMDASNESSetJacobianLocal");
99
100 ierr = DMSetMatType(*dm, "baij");
101 PISM_CHK(ierr, "DMSetMatType");
102
103 ierr = DMSetApplicationContext(*dm, &m_callback_data);
104 PISM_CHK(ierr, "DMSetApplicationContext");
105
106 ierr = SNESSetDM(m_snes, *dm);
107 PISM_CHK(ierr, "SNESSetDM");
108
109 // Default of maximum 200 iterations; possibly overridden by command line options
110 int snes_max_it = 200;
111 ierr = SNESSetTolerances(m_snes, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, snes_max_it,
112 PETSC_DEFAULT);
113 PISM_CHK(ierr, "SNESSetTolerances");
114
115 ierr = SNESSetFromOptions(m_snes);
116 PISM_CHK(ierr, "SNESSetFromOptions");
117
119 "node types: interior, boundary, exterior"); // no units or standard name
120
121 // Element::nodal_values() expects a ghosted array::Scalar. Ghosts if this field are
122 // never assigned to and not communicated, though.
124 "residual contribution from lateral boundaries"); // no units or standard name
125}
126
127// Initialize the solver, called once by the client before use.
129
131
132 // Use explicit driving stress if provided.
133 if (m_grid->variables().is_available("ssa_driving_stress_x") and
134 m_grid->variables().is_available("ssa_driving_stress_y")) {
135 m_driving_stress_x = m_grid->variables().get_2d_scalar("ssa_driving_stress_x");
136 m_driving_stress_y = m_grid->variables().get_2d_scalar("ssa_driving_stress_y");
137 }
138
139 m_log->message(2, " [using the SNES-based finite element method implementation]\n");
140
141 // process command-line options
142 {
143 m_dirichletScale = 1.0e9;
144 m_dirichletScale = options::Real(m_sys, "-ssa_fe_dirichlet_scale",
145 "Enforce Dirichlet conditions with this additional scaling",
146 "1", m_dirichletScale);
147 }
148
149 // On restart, SSA::init() reads the SSA velocity from a PISM output file
150 // into array::Vector "velocity". We use that field as an initial guess.
151 // If we are not restarting from a PISM file, "velocity" is identically zero,
152 // and the call below clears m_velocity_global.
153
155}
156
157/** Solve the SSA system of equations.
158
159 The FEM solver computes values of the various coefficients (most notably: the vertically-averaged
160 ice hardness) at each of the element nodes *only once*.
161
162 When running in an ice-model context, at each time step, SSA::update() is called, which calls
163 SSAFEM::solve(). Since coefficients have generally changed between time steps, we need to recompute
164 coefficients. On the other hand, in the context of inversion, coefficients will not change between
165 iteration and there is no need to recompute their values.
166
167 So there are two different solve methods, SSAFEM::solve() and SSAFEM::solve_nocache(). The only
168 difference is that SSAFEM::solve() recomputes the cached values of the coefficients before calling
169 SSAFEM::solve_nocache().
170 */
171void SSAFEM::solve(const Inputs &inputs) {
172
173 auto reason = solve_with_reason(inputs);
174 if (reason->failed()) {
176 "SSAFEM solve failed to converge (SNES reason %s)",
177 reason->description().c_str());
178 }
179
180 if (m_log->get_threshold() > 2) {
181 m_stdout_ssa += "SSAFEM converged (SNES reason " + reason->description() + ")";
182 }
183}
184
185std::shared_ptr<TerminationReason> SSAFEM::solve_with_reason(const Inputs &inputs) {
186
187 // Set up the system to solve.
188 cache_inputs(inputs);
189
190 return solve_nocache();
191}
192
193//! Solve the SSA without first recomputing the values of coefficients at quad
194//! points. See the disccusion of SSAFEM::solve for more discussion.
195std::shared_ptr<TerminationReason> SSAFEM::solve_nocache() {
196 PetscErrorCode ierr;
197
198 m_epsilon_ssa = m_config->get_number("stress_balance.ssa.epsilon");
199
200 options::String filename("-ssa_view", "");
201 if (filename.is_set()) {
202 petsc::Viewer viewer;
203 ierr = PetscViewerASCIIOpen(m_grid->com, filename->c_str(), viewer.rawptr());
204 PISM_CHK(ierr, "PetscViewerASCIIOpen");
205
206 ierr = PetscViewerASCIIPrintf(viewer, "SNES before SSASolve_FE\n");
207 PISM_CHK(ierr, "PetscViewerASCIIPrintf");
208
209 ierr = SNESView(m_snes, viewer);
210 PISM_CHK(ierr, "SNESView");
211
212 ierr = PetscViewerASCIIPrintf(viewer, "solution vector before SSASolve_FE\n");
213 PISM_CHK(ierr, "PetscViewerASCIIPrintf");
214
215 ierr = VecView(m_velocity_global.vec(), viewer);
216 PISM_CHK(ierr, "VecView");
217 }
218
219 m_stdout_ssa.clear();
220 if (m_log->get_threshold() >= 2) {
221 m_stdout_ssa = " SSA: ";
222 }
223
224 // Solve:
225 ierr = SNESSolve(m_snes, NULL, m_velocity_global.vec());
226 PISM_CHK(ierr, "SNESSolve");
227
228 // See if it worked.
229 SNESConvergedReason snes_reason;
230 ierr = SNESGetConvergedReason(m_snes, &snes_reason); PISM_CHK(ierr, "SNESGetConvergedReason");
231
232 std::shared_ptr<TerminationReason> reason(new SNESTerminationReason(snes_reason));
233 if (not reason->failed()) {
234
235 // Extract the solution back from SSAX to velocity and communicate.
237
238 bool view_solution = options::Bool("-ssa_view_solution", "view solution of the SSA system");
239 if (view_solution) {
240 petsc::Viewer viewer;
241 ierr = PetscViewerASCIIOpen(m_grid->com, filename->c_str(), viewer.rawptr());
242 PISM_CHK(ierr, "PetscViewerASCIIOpen");
243
244 ierr = PetscViewerASCIIPrintf(viewer, "solution vector after SSASolve\n");
245 PISM_CHK(ierr, "PetscViewerASCIIPrintf");
246
247 ierr = VecView(m_velocity_global.vec(), viewer);
248 PISM_CHK(ierr, "VecView");
249 }
250
251 }
252
253 return reason;
254}
255
256//! Initialize stored data from the coefficients in the SSA. Called by SSAFEM::solve.
257/**
258 This method is should be called after SSAFEM::init and whenever any geometry or temperature
259 related coefficients have changed. The method stores the values of the coefficients the nodes of
260 each element so that these are available to the residual and Jacobian evaluation methods.
261
262 We store the vertical average of the ice hardness to avoid re-doing this computation every
263 iteration.
264
265 In addition to coefficients at element nodes we store "node types" used to identify interior
266 elements, exterior elements, and boundary faces.
267*/
268void SSAFEM::cache_inputs(const Inputs &inputs) {
269
270 // Make copies of BC mask and BC values: they are needed in SNES callbacks and
271 // inputs.bc_{mask,values} are not available there.
272 if ((inputs.bc_mask != nullptr) and (inputs.bc_values != nullptr)) {
273 m_bc_mask.copy_from(*inputs.bc_mask);
275 } else {
276 m_bc_mask.set(0.0);
277 }
278
279 const std::vector<double> &z = m_grid->z();
280
282 inputs.enthalpy,
283 &inputs.geometry->ice_thickness,
284 &inputs.geometry->bed_elevation,
286 inputs.basal_yield_stress};
287
288 bool use_explicit_driving_stress = (m_driving_stress_x != NULL) && (m_driving_stress_y != NULL);
289 if (use_explicit_driving_stress) {
291 }
292
293 ParallelSection loop(m_grid->com);
294 try {
295 for (auto p : m_grid->points()) {
296 const int i = p.i(), j = p.j();
297
298 double thickness = inputs.geometry->ice_thickness(i, j);
299
300 Vector2d tau_d;
301 if (use_explicit_driving_stress) {
302 tau_d.u = (*m_driving_stress_x)(i, j);
303 tau_d.v = (*m_driving_stress_y)(i, j);
304 } else {
305 // tau_d above is set to zero by the Vector2d
306 // constructor, but is not used
307 }
308
309 const double *enthalpy = inputs.enthalpy->get_column(i, j);
310 double hardness = rheology::averaged_hardness(*m_flow_law, thickness,
311 m_grid->kBelowHeight(thickness),
312 z.data(), enthalpy);
313
314 Coefficients c;
315 c.thickness = thickness;
316 c.bed = inputs.geometry->bed_elevation(i, j);
317 c.sea_level = inputs.geometry->sea_level_elevation(i, j);
318 c.tauc = (*inputs.basal_yield_stress)(i, j);
319 c.hardness = hardness;
320 c.driving_stress = tau_d;
321
322 m_coefficients(i, j) = c;
323 } // loop over owned grid points
324 } catch (...) {
325 loop.failed();
326 }
327 loop.check();
328
329 m_coefficients.update_ghosts();
330
331 const bool use_cfbc = m_config->get_flag("stress_balance.calving_front_stress_bc");
332 if (use_cfbc) {
333 // Note: the call below uses ghosts of inputs.geometry->ice_thickness.
335 m_config->get_number("stress_balance.ice_free_thickness_standard"),
337 } else {
339 }
340
341 cache_residual_cfbc(inputs);
342
343}
344
345//! Compute quadrature point values of various coefficients given a quadrature `Q` and nodal values.
347 const Coefficients *x,
348 int *mask,
349 double *thickness,
350 double *tauc,
351 double *hardness) const {
352 // quadrature size
353 unsigned int n = E.n_pts();
354
355 for (unsigned int q = 0; q < n; q++) {
356 double
357 bed = 0.0,
358 sea_level = 0.0;
359
360 thickness[q] = 0.0;
361 tauc[q] = 0.0;
362 hardness[q] = 0.0;
363
364 for (int k = 0; k < E.n_chi(); k++) {
365 const fem::Germ &psi = E.chi(q, k);
366
367 thickness[q] += psi.val * x[k].thickness;
368 bed += psi.val * x[k].bed;
369 sea_level += psi.val * x[k].sea_level;
370 tauc[q] += psi.val * x[k].tauc;
371 hardness[q] += psi.val * x[k].hardness;
372 }
373
374 mask[q] = m_gc.mask(sea_level, bed, thickness[q]);
375 }
376}
377
378//! Compute gravitational driving stress at quadrature points.
379//! Uses explicitly-provided nodal values.
381 const Coefficients *x,
382 Vector2d *result) {
383 const unsigned int n = E.n_pts();
384
385 for (unsigned int q = 0; q < n; q++) {
386 result[q] = 0.0;
387
388 for (int k = 0; k < E.n_chi(); k++) {
389 const fem::Germ &psi = E.chi(q, k);
390 result[q] += psi.val * x[k].driving_stress;
391 }
392 }
393}
394
395/** Compute the the driving stress at quadrature points.
396 The surface elevation @f$ h @f$ is defined by
397 @f{equation}{
398 h =
399 \begin{cases}
400 b + H, & \mathrm{grounded},\\
401 z_{sl} + \alpha H, & \mathrm{shelf},
402 \end{cases}
403 @f}
404 where @f$ b @f$ is the bed elevation, @f$ H @f$ is the ice thickness, and @f$ z_{sl} @f$ is the
405 sea level, and @f$ \alpha @f$ is the ratio of densities of ice and sea water.
406
407 Note that if @f$ b @f$, @f$ H @f$, and @f$ z_{sl} @f$ are bilinear (or @f$ C^{1} @f$ in
408 general), @f$ h @f$ is continuous but not continuously differentiable. In
409 other words, the gradient of @f$ h @f$ is not continuous, so near the
410 grounding line we cannot compute the gravitational driving stress
411 @f$ \tau_{d} = - \rho g H \nabla h @f$ using the @f$ Q_1 @f$ or @f$ P_1 @f$ FE basis
412 expansion of @f$ h @f$.
413
414 We differentiate the equation above instead, getting
415 @f{equation}{
416 \nabla h =
417 \begin{cases}
418 \nabla b + \nabla H, & \mathrm{grounded},\\
419 \nabla z_{sl} + \alpha \nabla H, & \mathrm{shelf}.
420 \end{cases}
421 @f}
422
423 and use basis expansions of @f$ \nabla b @f$ and @f$ \nabla H @f$.
424
425 Overall, we have
426
427 @f{align*}{
428 \tau_{d} &= \rho g H \nabla h\\
429 &=
430 - \rho g H
431 \begin{cases}
432 \nabla b + \nabla H, & \mathrm{grounded},\\
433 \alpha \nabla H, & \mathrm{shelf},
434 \end{cases}
435 @f}
436
437 because @f$ z = z_{sl} @f$ defines the geoid surface and so *its gradient
438 does not contribute to the driving stress*.
439*/
441 const Coefficients *x,
442 Vector2d *result) const {
443 const unsigned int n = E.n_pts();
444
445 for (unsigned int q = 0; q < n; q++) {
446 double
447 sea_level = 0.0,
448 b = 0.0,
449 b_x = 0.0,
450 b_y = 0.0,
451 H = 0.0,
452 H_x = 0.0,
453 H_y = 0.0;
454
455 result[q] = 0.0;
456
457 for (int k = 0; k < E.n_chi(); k++) {
458 const fem::Germ &psi = E.chi(q, k);
459
460 b += psi.val * x[k].bed;
461 b_x += psi.dx * x[k].bed;
462 b_y += psi.dy * x[k].bed;
463
464 H += psi.val * x[k].thickness;
465 H_x += psi.dx * x[k].thickness;
466 H_y += psi.dy * x[k].thickness;
467
468 sea_level += psi.val * x[k].sea_level;
469 }
470
471 const int M = m_gc.mask(sea_level, b, H);
472 const bool grounded = mask::grounded(M);
473
474 const double
475 pressure = m_rho_g * H,
476 h_x = grounded ? b_x + H_x : m_alpha * H_x,
477 h_y = grounded ? b_y + H_y : m_alpha * H_y;
478
479 result[q].u = - pressure * h_x;
480 result[q].v = - pressure * h_y;
481 }
482}
483
484
485/** @brief Compute the "(regularized effective viscosity) x (ice thickness)" and effective viscous
486 * bed strength from the current solution, at a single quadrature point.
487 *
488 * @param[in] thickness ice thickness
489 * @param[in] hardness ice hardness
490 * @param[in] mask cell type mask
491 * @param[in] tauc basal yield stress
492 * @param[in] U the value of the solution
493 * @param[in] U_x x-derivatives of velocity components
494 * @param[in] U_y y-derivatives of velocity components
495 * @param[out] nuH product of the ice viscosity and thickness @f$ \nu H @f$
496 * @param[out] dnuH derivative of @f$ \nu H @f$ with respect to the
497 * second invariant @f$ \gamma @f$. Set to NULL if
498 * not desired.
499 * @param[out] beta basal drag coefficient @f$ \beta @f$
500 * @param[out] dbeta derivative of @f$ \beta @f$ with respect to the
501 * second invariant @f$ \gamma @f$. Set to NULL if
502 * not desired.
503 */
504void SSAFEM::PointwiseNuHAndBeta(double thickness,
505 double hardness,
506 int mask,
507 double tauc,
508 const Vector2d &U,
509 const Vector2d &U_x,
510 const Vector2d &U_y,
511 double *nuH, double *dnuH,
512 double *beta, double *dbeta) {
513
514 if (thickness < strength_extension->get_min_thickness()) {
516 if (dnuH != nullptr) {
517 *dnuH = 0;
518 }
519 } else {
520 m_flow_law->effective_viscosity(hardness, secondInvariant_2D(U_x, U_y),
521 nuH, dnuH);
522
523 *nuH = m_epsilon_ssa + *nuH * thickness;
524
525 if (dnuH != nullptr) {
526 *dnuH *= thickness;
527 }
528 }
529
530 if (mask::grounded_ice(mask)) {
531 m_basal_sliding_law->drag_with_derivative(tauc, U.u, U.v, beta, dbeta);
532 } else {
533 *beta = 0;
534
535 if (mask::ice_free_land(mask)) {
537 }
538
539 if (dbeta != nullptr) {
540 *dbeta = 0;
541 }
542 }
543}
544
545
546//! Compute and cache residual contributions from the integral over the lateral boundary.
547/**
548
549 This method computes FIXME.
550
551 */
553
554 // Compute 1D (linear) shape function values at quadrature points on element boundaries.
555 // These values do not depend on the element type or dimensions of a particular physical
556 // element.
557 //
558 // Note: the number of quadrature points (2) is hard-wired below.
559 const unsigned int Nq = 2;
560 double chi_b[Nq][Nq];
561 {
562 // interval length does not matter here
563 fem::Gaussian2 Q(1.0);
564
565 for (int i = 0; i < 2; ++i) { // 2 functions
566 for (int j = 0; j < 2; ++j) { // 2 quadrature points
567 chi_b[i][j] = fem::linear::chi(i, Q.point(j)).val;
568 }
569 }
570 }
571
572 const unsigned int Nk = fem::q1::n_chi;
573
574 using fem::P1Element2;
576 P1Element2 p1_element[Nk] = {P1Element2(*m_grid, Q_p1, 0),
577 P1Element2(*m_grid, Q_p1, 1),
578 P1Element2(*m_grid, Q_p1, 2),
579 P1Element2(*m_grid, Q_p1, 3)};
580
581 using mask::ocean;
582
583 const bool
584 use_cfbc = m_config->get_flag("stress_balance.calving_front_stress_bc");
585
586 const double
587 ice_density = m_config->get_number("constants.ice.density"),
588 ocean_density = m_config->get_number("constants.sea_water.density"),
589 standard_gravity = m_config->get_number("constants.standard_gravity");
590
591 // Reset the boundary integral so that all values are overwritten.
593
594 if (not use_cfbc) {
595 // If CFBC is not used then we're done.
596 return;
597 }
598
600 &inputs.geometry->ice_thickness,
601 &inputs.geometry->bed_elevation,
604
605 // Iterate over the elements.
606 const int
607 xs = m_element_index.xs,
608 xm = m_element_index.xm,
609 ys = m_element_index.ys,
610 ym = m_element_index.ym;
611
612 ParallelSection loop(m_grid->com);
613 try {
614 for (int j = ys; j < ys + ym; j++) {
615 for (int i = xs; i < xs + xm; i++) {
616 // Initialize the map from global to element degrees of freedom.
617 m_q1_element.reset(i, j);
618
619 int node_type[Nk];
621
622 fem::Element2 *E = nullptr;
623 {
624 auto type = fem::element_type(node_type);
625
626 switch (type) {
628 continue;
629 case fem::ELEMENT_Q:
630 E = &m_q1_element;
631 break;
632 default:
633 E = &p1_element[type];
634
635 E->reset(i, j);
636
637 // get node types *again*, this time at the nodes of this P1 element
638 E->nodal_values(m_node_type, node_type);
639 }
640 }
641
642 // residual contributions at element nodes
643 std::vector<Vector2d> I(Nk);
644
645 double H_nodal[Nk];
646 E->nodal_values(inputs.geometry->ice_thickness.array(), H_nodal);
647
648 double b_nodal[Nk];
649 E->nodal_values(inputs.geometry->bed_elevation.array(), b_nodal);
650
651 double sl_nodal[Nk];
652 E->nodal_values(inputs.geometry->sea_level_elevation.array(), sl_nodal);
653
654 // storage for values of test functions on a side of the element
655 double psi[2] = {0.0, 0.0};
656
657 const unsigned int n_sides = E->n_sides();
658 // loop over element sides
659 for (unsigned int s = 0; s < n_sides; ++s) {
660
661 // initialize the quadrature and compute quadrature weights
663
664 // nodes incident to the current side
665 const int
666 n0 = s,
667 n1 = (n0 + 1) % n_sides;
668
669 if (not (node_type[n0] == NODE_BOUNDARY and
670 node_type[n1] == NODE_BOUNDARY)) {
671 // not a boundary side; skip it
672 continue;
673 }
674
675 for (unsigned int q = 0; q < Nq; ++q) {
676
677 double W = Q.weight(q);
678
679 // test functions at nodes incident to the current side, evaluated at the
680 // quadrature point q
681 psi[0] = chi_b[0][q];
682 psi[1] = chi_b[1][q];
683
684 // Compute ice thickness and bed elevation at a quadrature point. This uses a 1D basis
685 // expansion on the current side.
686 const double
687 H = H_nodal[n0] * psi[0] + H_nodal[n1] * psi[1],
688 bed = b_nodal[n0] * psi[0] + b_nodal[n1] * psi[1],
689 sea_level = sl_nodal[n0] * psi[0] + sl_nodal[n1] * psi[1];
690
691 // vertically averaged pressures at a quadrature point
692 double
693 P_ice = 0.5 * ice_density * standard_gravity * H,
694 P_water = average_water_column_pressure(H, bed, sea_level,
695 ice_density, ocean_density,
696 standard_gravity);
697
698 // vertical integral of the pressure difference
699 double dP = H * (P_ice - P_water);
700 // FIXME: implement melange pressure forcing
701
702 // This integral contributes to the residual at 2 nodes (the ones incident to the
703 // current side). This is is written in a way that allows *adding* (... += ...) the
704 // boundary contribution in the residual computation.
705 //
706 // FIXME: I need to include the special case corresponding to ice margins next
707 // to fjord walls, nunataks, etc. In this case dP == 0.
708 //
709 // FIXME: set pressure difference to zero at grounded locations at domain
710 // boundaries.
711 I[n0] += W * (- psi[0] * dP) * E->normal(s);
712 I[n1] += W * (- psi[1] * dP) * E->normal(s);
713 } // q-loop
714
715 } // loop over element sides
716
718
719 } // i-loop
720 } // j-loop
721 } catch (...) {
722 loop.failed();
723 }
724 loop.check();
725}
726
727
728//! Implements the callback for computing the residual.
729/*!
730 * Compute the residual \f[r_{ij}= G(x, \psi_{ij}) \f] where \f$G\f$
731 * is the weak form of the SSA, \f$x\f$ is the current approximate
732 * solution, and the \f$\psi_{ij}\f$ are test functions.
733 *
734 * The weak form of the SSA system is
735 */
736void SSAFEM::compute_local_function(Vector2d const *const *const velocity_global,
737 Vector2d **residual_global) {
738
739 const bool use_explicit_driving_stress = (m_driving_stress_x != NULL) && (m_driving_stress_y != NULL);
740
741 const bool use_cfbc = m_config->get_flag("stress_balance.calving_front_stress_bc");
742
743 const unsigned int Nk = fem::q1::n_chi;
744 const unsigned int Nq_max = fem::MAX_QUADRATURE_SIZE;
745
746 using fem::P1Element2;
748 P1Element2 p1_element[Nk] = {P1Element2(*m_grid, Q_p1, 0),
749 P1Element2(*m_grid, Q_p1, 1),
750 P1Element2(*m_grid, Q_p1, 2),
751 P1Element2(*m_grid, Q_p1, 3)};
752
754
755 // Set the boundary contribution of the residual. This is computed at the nodes, so we don't want
756 // to set it using Element::add_contribution() because that would lead to
757 // double-counting. Also note that without CFBC m_boundary_integral is exactly zero.
758 for (auto p : m_grid->points()) {
759 const int i = p.i(), j = p.j();
760
761 residual_global[j][i] = m_boundary_integral(i, j);
762 }
763
764 // Start access to Dirichlet data if present.
766
767 // Storage for the current solution and its derivatives at quadrature points.
768 Vector2d U[Nq_max], U_x[Nq_max], U_y[Nq_max];
769
770 // Iterate over the elements.
771 const int
772 xs = m_element_index.xs,
773 xm = m_element_index.xm,
774 ys = m_element_index.ys,
775 ym = m_element_index.ym;
776
777 ParallelSection loop(m_grid->com);
778 try {
779 for (int j = ys; j < ys + ym; j++) {
780 for (int i = xs; i < xs + xm; i++) {
781
782 fem::Element2 *E = nullptr;
783 {
784 m_q1_element.reset(i, j);
785
786 int node_type[Nk];
788
789 auto type = fem::element_type(node_type);
790
791 if (use_cfbc) {
792 if (type == fem::ELEMENT_EXTERIOR) {
793 // skip exterior elements
794 continue;
795 }
796
797 if (type == fem::ELEMENT_Q) {
798 E = &m_q1_element;
799 } else {
800 E = &p1_element[type];
801
802 E->reset(i, j);
803 }
804 } else {
805 // if use_cfbc == false all elements are interior and Q1
806 E = &m_q1_element;
807 }
808 }
809
810 // Number of quadrature points.
811 const unsigned int Nq = E->n_pts();
812
813 // Storage for the solution and residuals at element nodes.
814 Vector2d residual[Nk];
815
816 int mask[Nq_max];
817 double thickness[Nq_max];
818 double tauc[Nq_max];
819 double hardness[Nq_max];
820 Vector2d tau_d[Nq_max];
821
822 {
823 Coefficients coeffs[Nk];
824 E->nodal_values(m_coefficients.array(), coeffs);
825
826 quad_point_values(*E, coeffs, mask, thickness, tauc, hardness);
827
828 if (use_explicit_driving_stress) {
829 explicit_driving_stress(*E, coeffs, tau_d);
830 } else {
831 driving_stress(*E, coeffs, tau_d);
832 }
833 }
834
835 {
836 // Obtain the value of the solution at the nodes
837 Vector2d velocity_nodal[Nk];
838 E->nodal_values(velocity_global, velocity_nodal);
839
840 // These values now need to be adjusted if some nodes in the element have Dirichlet data.
841 if (dirichlet_data) {
842 // Set elements of velocity_nodal that correspond to Dirichlet nodes to prescribed
843 // values.
844 dirichlet_data.enforce(*E, velocity_nodal);
845 // mark Dirichlet nodes in E so that they are not touched by
846 // add_contribution() below
847 dirichlet_data.constrain(*E);
848 }
849
850 // Compute the solution values and its gradient at the quadrature points.
851 E->evaluate(velocity_nodal, // input
852 U, U_x, U_y); // outputs
853 }
854
855 // Zero out the element-local residual in preparation for updating it.
856 for (unsigned int k = 0; k < Nk; k++) {
857 residual[k] = 0;
858 }
859
860 // loop over quadrature points:
861 for (unsigned int q = 0; q < Nq; q++) {
862
863 auto W = E->weight(q);
864
865 double eta = 0.0, beta = 0.0;
866 PointwiseNuHAndBeta(thickness[q], hardness[q], mask[q], tauc[q],
867 U[q], U_x[q], U_y[q], // inputs
868 &eta, NULL, &beta, NULL); // outputs
869
870 // The next few lines compute the actual residual for the element.
871 const Vector2d tau_b = U[q] * (- beta); // basal shear stress
872
873 const double
874 u_x = U_x[q].u,
875 v_y = U_y[q].v,
876 u_y_plus_v_x = U_y[q].u + U_x[q].v;
877
878 // Loop over test functions.
879 for (int k = 0; k < E->n_chi(); k++) {
880 const fem::Germ &psi = E->chi(q, k);
881
882 residual[k].u += W * (eta * (psi.dx * (4.0 * u_x + 2.0 * v_y) + psi.dy * u_y_plus_v_x)
883 - psi.val * (tau_b.u + tau_d[q].u));
884 residual[k].v += W * (eta * (psi.dx * u_y_plus_v_x + psi.dy * (2.0 * u_x + 4.0 * v_y))
885 - psi.val * (tau_b.v + tau_d[q].v));
886 } // k (test functions)
887 } // q (quadrature points)
888
889 E->add_contribution(residual, residual_global);
890 } // i-loop
891 } // j-loop
892 } catch (...) {
893 loop.failed();
894 }
895 loop.check();
896
897 // Until now we have not touched rows in the residual corresponding to Dirichlet data.
898 // We fix this now.
899 if (dirichlet_data) {
900 dirichlet_data.fix_residual(velocity_global, residual_global);
901 }
902
903 if (use_cfbc) {
904 // Prescribe homogeneous Dirichlet B.C. at ice-free nodes. This uses the fact that m_node_type
905 // can be used as a "Dirichlet B.C. mask", i.e. ice-free nodes (and only ice-free nodes) are
906 // marked with ones.
907 fem::DirichletData_Vector dirichlet_ice_free(&m_node_type, NULL, m_dirichletScale);
908 dirichlet_ice_free.fix_residual_homogeneous(residual_global);
909 }
910
911 monitor_function(velocity_global, residual_global);
912}
913
914void SSAFEM::monitor_function(Vector2d const *const *const velocity_global,
915 Vector2d const *const *const residual_global) {
916 PetscErrorCode ierr;
917 bool monitorFunction = options::Bool("-ssa_monitor_function", "monitor the SSA residual");
918 if (not monitorFunction) {
919 return;
920 }
921
922 ierr = PetscPrintf(m_grid->com,
923 "SSA Solution and Function values (pointwise residuals)\n");
924 PISM_CHK(ierr, "PetscPrintf");
925
926 ParallelSection loop(m_grid->com);
927 try {
928 for (auto p : m_grid->points()) {
929 const int i = p.i(), j = p.j();
930
931 ierr = PetscSynchronizedPrintf(m_grid->com,
932 "[%2d, %2d] u=(%12.10e, %12.10e) f=(%12.4e, %12.4e)\n",
933 i, j,
934 velocity_global[j][i].u, velocity_global[j][i].v,
935 residual_global[j][i].u, residual_global[j][i].v);
936 PISM_CHK(ierr, "PetscSynchronizedPrintf");
937 }
938 } catch (...) {
939 loop.failed();
940 }
941 loop.check();
942
943 ierr = PetscSynchronizedFlush(m_grid->com, NULL);
944 PISM_CHK(ierr, "PetscSynchronizedFlush");
945}
946
947
948//! Implements the callback for computing the Jacobian.
949/*!
950 Compute the Jacobian
951
952 @f[ J_{ij}{kl} \frac{d r_{ij}}{d x_{kl}}= G(x, \psi_{ij}) @f]
953
954 where \f$G\f$ is the weak form of the SSA, \f$x\f$ is the current
955 approximate solution, and the \f$\psi_{ij}\f$ are test functions.
956
957*/
958void SSAFEM::compute_local_jacobian(Vector2d const *const *const velocity_global, Mat Jac) {
959
960 const bool use_cfbc = m_config->get_flag("stress_balance.calving_front_stress_bc");
961
962 const unsigned int Nk = fem::q1::n_chi;
963 const unsigned int Nq_max = fem::MAX_QUADRATURE_SIZE;
964
965 using fem::P1Element2;
967 P1Element2 p1_element[Nk] = {P1Element2(*m_grid, Q_p1, 0),
968 P1Element2(*m_grid, Q_p1, 1),
969 P1Element2(*m_grid, Q_p1, 2),
970 P1Element2(*m_grid, Q_p1, 3)};
971
972 // Zero out the Jacobian in preparation for updating it.
973 PetscErrorCode ierr = MatZeroEntries(Jac);
974 PISM_CHK(ierr, "MatZeroEntries");
975
977
978 // Start access to Dirichlet data if present.
980
981 // Storage for the current solution at quadrature points.
982 Vector2d U[Nq_max], U_x[Nq_max], U_y[Nq_max];
983
984 // Loop through all the elements.
985 int
986 xs = m_element_index.xs,
987 xm = m_element_index.xm,
988 ys = m_element_index.ys,
989 ym = m_element_index.ym;
990
991 ParallelSection loop(m_grid->com);
992 try {
993 for (int j = ys; j < ys + ym; j++) {
994 for (int i = xs; i < xs + xm; i++) {
995
996 fem::Element2 *E = nullptr;
997 {
998 m_q1_element.reset(i, j);
999
1000 int node_type[Nk];
1002
1003 auto type = fem::element_type(node_type);
1004
1005 if (use_cfbc) {
1006 if (type == fem::ELEMENT_EXTERIOR) {
1007 // skip exterior elements
1008 continue;
1009 }
1010
1011 if (type == fem::ELEMENT_Q) {
1012 E = &m_q1_element;
1013 } else {
1014 E = &p1_element[type];
1015
1016 E->reset(i, j);
1017 }
1018 } else {
1019 // if use_cfbc == false all elements are interior and Q1
1020 E = &m_q1_element;
1021 }
1022 }
1023
1024 // Number of quadrature points.
1025 const unsigned int
1026 Nq = E->n_pts(),
1027 n_chi = E->n_chi();
1028
1029 int mask[Nq_max];
1030 double thickness[Nq_max];
1031 double tauc[Nq_max];
1032 double hardness[Nq_max];
1033
1034 {
1035 Coefficients coeffs[Nk];
1036 E->nodal_values(m_coefficients.array(), coeffs);
1037
1038 quad_point_values(*E, coeffs,
1039 mask, thickness, tauc, hardness);
1040 }
1041
1042 {
1043 // Values of the solution at the nodes of the current element.
1044 Vector2d velocity_nodal[Nk];
1045 E->nodal_values(velocity_global, velocity_nodal);
1046
1047 // These values now need to be adjusted if some nodes in the element have
1048 // Dirichlet data.
1049 if (dirichlet_data) {
1050 dirichlet_data.enforce(*E, velocity_nodal);
1051 dirichlet_data.constrain(*E);
1052 }
1053 // Compute the values of the solution at the quadrature points.
1054 E->evaluate(velocity_nodal, U, U_x, U_y);
1055 }
1056
1057 // Element-local Jacobian matrix (there are Nk vector valued degrees
1058 // of freedom per element, for a total of (2*Nk)*(2*Nk) = 64
1059 // entries in the local Jacobian.
1060 double K[2*Nk][2*Nk];
1061 // Build the element-local Jacobian.
1062 ierr = PetscMemzero(K, sizeof(K));
1063 PISM_CHK(ierr, "PetscMemzero");
1064
1065 for (unsigned int q = 0; q < Nq; q++) {
1066
1067 const double
1068 W = E->weight(q),
1069 u = U[q].u,
1070 v = U[q].v,
1071 u_x = U_x[q].u,
1072 v_y = U_y[q].v,
1073 u_y_plus_v_x = U_y[q].u + U_x[q].v;
1074
1075 double eta = 0.0, deta = 0.0, beta = 0.0, dbeta = 0.0;
1076 PointwiseNuHAndBeta(thickness[q], hardness[q], mask[q], tauc[q],
1077 U[q], U_x[q], U_y[q],
1078 &eta, &deta, &beta, &dbeta);
1079
1080 for (unsigned int l = 0; l < n_chi; l++) { // Trial functions
1081
1082 // Current trial function and its derivatives:
1083 const fem::Germ &phi = E->chi(q, l);
1084
1085 // Derivatives of \gamma with respect to u_l and v_l:
1086 const double
1087 gamma_u = (2.0 * u_x + v_y) * phi.dx + 0.5 * u_y_plus_v_x * phi.dy,
1088 gamma_v = 0.5 * u_y_plus_v_x * phi.dx + (u_x + 2.0 * v_y) * phi.dy;
1089
1090 // Derivatives of \eta = \nu*H with respect to u_l and v_l:
1091 const double
1092 eta_u = deta * gamma_u,
1093 eta_v = deta * gamma_v;
1094
1095 // Derivatives of the basal shear stress term (\tau_b):
1096 const double
1097 taub_xu = -dbeta * u * u * phi.val - beta * phi.val, // x-component, derivative with respect to u_l
1098 taub_xv = -dbeta * u * v * phi.val, // x-component, derivative with respect to u_l
1099 taub_yu = -dbeta * v * u * phi.val, // y-component, derivative with respect to v_l
1100 taub_yv = -dbeta * v * v * phi.val - beta * phi.val; // y-component, derivative with respect to v_l
1101
1102 for (unsigned int k = 0; k < n_chi; k++) { // Test functions
1103
1104 // Current test function and its derivatives:
1105
1106 const fem::Germ &psi = E->chi(q, k);
1107
1108 if (eta == 0) {
1109 ierr = PetscPrintf(PETSC_COMM_SELF, "eta=0 i %d j %d q %d k %d\n", i, j, q, k);
1110 PISM_CHK(ierr, "PetscPrintf");
1111 }
1112
1113 // u-u coupling
1114 K[k*2 + 0][l*2 + 0] += W * (eta_u * (psi.dx * (4 * u_x + 2 * v_y) + psi.dy * u_y_plus_v_x)
1115 + eta * (4 * psi.dx * phi.dx + psi.dy * phi.dy) - psi.val * taub_xu);
1116 // u-v coupling
1117 K[k*2 + 0][l*2 + 1] += W * (eta_v * (psi.dx * (4 * u_x + 2 * v_y) + psi.dy * u_y_plus_v_x)
1118 + eta * (2 * psi.dx * phi.dy + psi.dy * phi.dx) - psi.val * taub_xv);
1119 // v-u coupling
1120 K[k*2 + 1][l*2 + 0] += W * (eta_u * (psi.dx * u_y_plus_v_x + psi.dy * (2 * u_x + 4 * v_y))
1121 + eta * (psi.dx * phi.dy + 2 * psi.dy * phi.dx) - psi.val * taub_yu);
1122 // v-v coupling
1123 K[k*2 + 1][l*2 + 1] += W * (eta_v * (psi.dx * u_y_plus_v_x + psi.dy * (2 * u_x + 4 * v_y))
1124 + eta * (psi.dx * phi.dx + 4 * psi.dy * phi.dy) - psi.val * taub_yv);
1125
1126 } // l
1127 } // k
1128 } // q
1129 E->add_contribution(&K[0][0], Jac);
1130 } // j
1131 } // i
1132 } catch (...) {
1133 loop.failed();
1134 }
1135 loop.check();
1136
1137
1138 // Until now, the rows and columns correspoinding to Dirichlet data
1139 // have not been set. We now put an identity block in for these
1140 // unknowns. Note that because we have takes steps to not touching
1141 // these columns previously, the symmetry of the Jacobian matrix is
1142 // preserved.
1143 if (dirichlet_data) {
1144 dirichlet_data.fix_jacobian(Jac);
1145 }
1146
1147 if (use_cfbc) {
1148 // Prescribe homogeneous Dirichlet B.C. at ice-free nodes. This uses the fact that m_node_type
1149 // can be used as a "Dirichlet B.C. mask", i.e. ice-free nodes (and only ice-free nodes) are
1150 // marked with ones.
1151 fem::DirichletData_Vector dirichlet_ice_free(&m_node_type, NULL, m_dirichletScale);
1152 dirichlet_ice_free.fix_jacobian(Jac);
1153 }
1154
1155 ierr = MatAssemblyBegin(Jac, MAT_FINAL_ASSEMBLY);
1156 PISM_CHK(ierr, "MatAssemblyBegin");
1157
1158 ierr = MatAssemblyEnd(Jac, MAT_FINAL_ASSEMBLY);
1159 PISM_CHK(ierr, "MatAssemblyEnd");
1160
1161 ierr = MatSetOption(Jac, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE);
1162 PISM_CHK(ierr, "MatSetOption");
1163
1164 ierr = MatSetOption(Jac, MAT_SYMMETRIC, PETSC_TRUE);
1165 PISM_CHK(ierr, "MatSetOption");
1166
1167 monitor_jacobian(Jac);
1168}
1169
1171 PetscErrorCode ierr;
1172 bool mon_jac = options::Bool("-ssa_monitor_jacobian", "monitor the SSA Jacobian");
1173
1174 if (not mon_jac) {
1175 return;
1176 }
1177
1178 // iter has to be a PetscInt because it is used in the
1179 // SNESGetIterationNumber() call below.
1180 PetscInt iter = 0;
1181 ierr = SNESGetIterationNumber(m_snes, &iter);
1182 PISM_CHK(ierr, "SNESGetIterationNumber");
1183
1184 auto file_name = pism::printf("PISM_SSAFEM_J%d.m", (int)iter);
1185
1186 m_log->message(2,
1187 "writing Matlab-readable file for SSAFEM system A xsoln = rhs to file `%s' ...\n",
1188 file_name.c_str());
1189
1190 petsc::Viewer viewer(m_grid->com);
1191
1192 ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);
1193 PISM_CHK(ierr, "PetscViewerSetType");
1194
1195 ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
1196 PISM_CHK(ierr, "PetscViewerPushFormat");
1197
1198 ierr = PetscViewerFileSetName(viewer, file_name.c_str());
1199 PISM_CHK(ierr, "PetscViewerFileSetName");
1200
1201 ierr = PetscObjectSetName((PetscObject) Jac, "A");
1202 PISM_CHK(ierr, "PetscObjectSetName");
1203
1204 ierr = MatView(Jac, viewer);
1205 PISM_CHK(ierr, "MatView");
1206
1207 ierr = PetscViewerPopFormat(viewer);
1208 PISM_CHK(ierr, "PetscViewerPopFormat");
1209}
1210
1211//!
1212PetscErrorCode SSAFEM::function_callback(DMDALocalInfo *info,
1213 Vector2d const *const *const velocity,
1214 Vector2d **residual,
1215 CallbackData *fe) {
1216 try {
1217 (void) info;
1218 fe->ssa->compute_local_function(velocity, residual);
1219 } catch (...) {
1220 MPI_Comm com = MPI_COMM_SELF;
1221 PetscErrorCode ierr = PetscObjectGetComm((PetscObject)fe->da, &com); CHKERRQ(ierr);
1223 SETERRQ(com, 1, "A PISM callback failed");
1224 }
1225 return 0;
1226}
1227
1228PetscErrorCode SSAFEM::jacobian_callback(DMDALocalInfo *info,
1229 Vector2d const *const *const velocity,
1230 Mat A, Mat J, CallbackData *fe) {
1231 try {
1232 (void) A;
1233 (void) info;
1235 } catch (...) {
1236 MPI_Comm com = MPI_COMM_SELF;
1237 PetscErrorCode ierr = PetscObjectGetComm((PetscObject)fe->da, &com); CHKERRQ(ierr);
1239 SETERRQ(com, 1, "A PISM callback failed");
1240 }
1241 return 0;
1242}
1243
1244} // end of namespace stressbalance
1245} // 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
std::shared_ptr< const Logger > m_log
logger (for easy access)
Definition Component.hh:164
int mask(double sea_level, double bed, double thickness) const
Definition Mask.hh:138
array::Scalar1 sea_level_elevation
Definition Geometry.hh:48
array::Scalar2 ice_thickness
Definition Geometry.hh:51
array::Scalar2 bed_elevation
Definition Geometry.hh:47
virtual void drag_with_derivative(double tauc, double vx, double vy, double *drag, double *ddrag) const
Compute the drag coefficient and its derivative with respect to .
void failed()
Indicates a failure of a parallel section.
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)
This class represents a 2D vector field (such as ice velocity) at a certain grid point.
Definition Vector2d.hh:29
T * rawptr()
Definition Wrapper.hh:34
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
void add(double alpha, const Array2D< T > &x)
Definition Array2D.hh:93
double * get_column(int i, int j)
Definition Array3D.cc:127
void set_interpolation_type(InterpolationType type)
Definition Array.cc:181
petsc::Vec & vec() const
Definition Array.cc:313
void set(double c)
Result: v[j] <- c for all j.
Definition Array.cc:659
std::shared_ptr< petsc::DM > dm() const
Definition Array.cc:327
VariableMetadata & metadata(unsigned int N=0)
Returns a reference to the VariableMetadata object containing metadata for the compoment N.
Definition Array.cc:467
void fix_residual(Vector2d const *const *x_global, Vector2d **r_global)
void enforce(const Element2 &element, Vector2d *x_e)
void fix_residual_homogeneous(Vector2d **r)
void constrain(Element2 &element)
Constrain element, i.e. ensure that quadratures do not contribute to Dirichlet nodes by marking corre...
void evaluate(const T *x, T *vals, T *dx, T *dy)
Given nodal values, compute the values and partial derivatives at the quadrature points.
Definition Element.hh:195
void reset(int i, int j)
Initialize the Element to element (i, j) for the purposes of inserting into global residual and Jacob...
Definition Element.cc:196
double side_length(unsigned int side) const
Definition Element.hh:185
void add_contribution(const T *local, T **y_global) const
Add the values of element-local contributions y to the global vector y_global.
Definition Element.hh:238
void nodal_values(const array::Scalar &x_global, int *result) const
Get nodal values of an integer mask.
Definition Element.cc:185
unsigned int n_sides() const
Definition Element.hh:181
Vector2d normal(unsigned int side) const
Definition Element.hh:176
int xm
total number of elements to loop over in the x-direction.
int ym
total number of elements to loop over in the y-direction.
int ys
y-coordinate of the first element to loop over.
int xs
x-coordinate of the first element to loop over.
int n_chi() const
Definition Element.hh:65
double weight(unsigned int q) const
Weight of the quadrature point q
Definition Element.hh:85
const Germ & chi(unsigned int q, unsigned int k) const
Definition Element.hh:73
unsigned int n_pts() const
Number of quadrature points.
Definition Element.hh:80
The mapping from global to local degrees of freedom.
Definition Element.hh:61
P1 element embedded in Q1Element2.
Definition Element.hh:273
double weight(int k) const
Definition Quadrature.hh:72
QuadPoint point(int k) const
Definition Quadrature.hh:68
bool is_set() const
Definition options.hh:35
const array::Scalar * basal_yield_stress
const array::Scalar * bc_mask
const array::Array3D * enthalpy
const array::Vector * bc_values
std::shared_ptr< TerminationReason > solve_with_reason(const Inputs &inputs)
Definition SSAFEM.cc:185
void PointwiseNuHAndBeta(double thickness, double hardness, int mask, double tauc, const Vector2d &U, const Vector2d &U_x, const Vector2d &U_y, double *nuH, double *dnuH, double *beta, double *dbeta)
Compute the "(regularized effective viscosity) x (ice thickness)" and effective viscous bed strength ...
Definition SSAFEM.cc:504
void driving_stress(const fem::Element &E, const Coefficients *x, Vector2d *result) const
Definition SSAFEM.cc:440
const array::Scalar * m_driving_stress_x
Definition SSAFEM.hh:142
void cache_inputs(const Inputs &inputs)
Initialize stored data from the coefficients in the SSA. Called by SSAFEM::solve.
Definition SSAFEM.cc:268
static PetscErrorCode function_callback(DMDALocalInfo *info, Vector2d const *const *velocity, Vector2d **residual, CallbackData *fe)
SNES callbacks.
Definition SSAFEM.cc:1212
SSAFEM(std::shared_ptr< const Grid > g)
Definition SSAFEM.cc:48
static void explicit_driving_stress(const fem::Element &E, const Coefficients *x, Vector2d *result)
Definition SSAFEM.cc:380
CallbackData m_callback_data
Definition SSAFEM.hh:122
void compute_local_jacobian(Vector2d const *const *velocity, Mat J)
Implements the callback for computing the Jacobian.
Definition SSAFEM.cc:958
array::Vector1 m_boundary_integral
Boundary integral (CFBC contribution to the residual).
Definition SSAFEM.hh:129
static PetscErrorCode jacobian_callback(DMDALocalInfo *info, Vector2d const *const *velocity, Mat A, Mat J, CallbackData *fe)
Definition SSAFEM.cc:1228
void monitor_jacobian(Mat Jac)
Definition SSAFEM.cc:1170
void monitor_function(Vector2d const *const *velocity_global, Vector2d const *const *residual_global)
Definition SSAFEM.cc:914
array::Scalar1 m_bc_mask
Definition SSAFEM.hh:67
virtual void solve(const Inputs &inputs)
Definition SSAFEM.cc:171
fem::Q1Element2 m_q1_element
Definition SSAFEM.hh:136
GeometryCalculator m_gc
Definition SSAFEM.hh:70
virtual void init_impl()
Initialize a generic regular-grid SSA solver.
Definition SSAFEM.cc:128
std::shared_ptr< TerminationReason > solve_nocache()
Definition SSAFEM.cc:195
void cache_residual_cfbc(const Inputs &inputs)
Compute and cache residual contributions from the integral over the lateral boundary.
Definition SSAFEM.cc:552
void compute_local_function(Vector2d const *const *velocity, Vector2d **residual)
Implements the callback for computing the residual.
Definition SSAFEM.cc:736
fem::ElementIterator m_element_index
Definition SSAFEM.hh:135
void quad_point_values(const fem::Element &E, const Coefficients *x, int *mask, double *thickness, double *tauc, double *hardness) const
Compute quadrature point values of various coefficients given a quadrature Q and nodal values.
Definition SSAFEM.cc:346
const array::Scalar * m_driving_stress_y
Definition SSAFEM.hh:143
array::Scalar1 m_node_type
Storage for node types (interior, boundary, exterior).
Definition SSAFEM.hh:127
array::Vector1 m_bc_values
Definition SSAFEM.hh:68
array::Array2D< Coefficients > m_coefficients
Definition SSAFEM.hh:74
double get_notional_strength() const
Returns strength = (viscosity times thickness).
Definition SSA.cc:61
array::Vector m_velocity_global
Definition SSA.hh:122
SSAStrengthExtension * strength_extension
Definition SSA.hh:104
std::string m_stdout_ssa
Definition SSA.hh:120
virtual void init_impl()
Initialize a generic regular-grid SSA solver.
Definition SSA.cc:104
PISM's SSA solver.
Definition SSA.hh:99
const array::Vector1 & velocity() const
Get the thickness-advective 2D velocity.
std::shared_ptr< rheology::FlowLaw > m_flow_law
IceBasalResistancePlasticLaw * m_basal_sliding_law
#define PISM_CHK(errcode, name)
#define PISM_ERROR_LOCATION
const double phi
Definition exactTestK.c:41
#define n
Definition exactTestM.c:37
Germ chi(unsigned int k, const QuadPoint &pt)
Linear basis functions on the interval [-1, -1].
Definition FEM.cc:30
const int n_chi
Definition FEM.hh:191
const unsigned int MAX_QUADRATURE_SIZE
Definition FEM.hh:173
ElementType element_type(const int node_type[q1::n_chi])
Definition FEM.cc:106
@ ELEMENT_Q
Definition FEM.hh:203
@ ELEMENT_EXTERIOR
Definition FEM.hh:205
bool grounded_ice(int M)
Definition Mask.hh:51
bool ice_free_land(int M)
Definition Mask.hh:64
bool grounded(int M)
Grounded cell (grounded ice or ice-free).
Definition Mask.hh:44
bool ocean(int M)
An ocean cell (floating ice or ice-free).
Definition Mask.hh:40
bool Bool(const std::string &option, const std::string &description)
Definition options.cc:190
double averaged_hardness(const FlowLaw &ice, double thickness, unsigned int kbelowH, const double *zlevels, const double *enthalpy)
Computes vertical average of B(E, p) ice hardness, namely .
Definition FlowLaw.cc:221
@ NODE_BOUNDARY
Definition node_types.hh:35
@ NODE_INTERIOR
Definition node_types.hh:34
static double secondInvariant_2D(const Vector2d &U_x, const Vector2d &U_y)
Definition FlowLaw.hh:45
void compute_node_types(const array::Scalar1 &ice_thickness, double thickness_threshold, array::Scalar &result)
Definition node_types.cc:62
double average_water_column_pressure(double ice_thickness, double bed, double floatation_level, double rho_ice, double rho_water, double g)
std::string printf(const char *format,...)
static const double k
Definition exactTestP.cc:42
void handle_fatal_errors(MPI_Comm com)
double dy
Function derivative with respect to y.
Definition FEM.hh:157
double val
Function value.
Definition FEM.hh:153
double dx
Function deriviative with respect to x.
Definition FEM.hh:155
Struct for gathering the value and derivative of a function at a point.
Definition FEM.hh:151
Adaptor for gluing SNESDAFormFunction callbacks to an SSAFEM.
Definition SSAFEM.hh:116
double tauc
basal yield stress
Definition SSAFEM.hh:60
Vector2d driving_stress
prescribed gravitational driving stress
Definition SSAFEM.hh:64