Loading [MathJax]/extensions/tex2jax.js
PISM, A Parallel Ice Sheet Model 2.2.2-d6b3a29ca committed by Constantine Khrulev on 2025-03-28
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
SSAFEM.hh
Go to the documentation of this file.
1// Copyright (C) 2009--2017, 2020, 2021, 2022, 2023, 2024, 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#ifndef _SSAFEM_H_
20#define _SSAFEM_H_
21
22#include "pism/stressbalance/ssa/SSA.hh"
23#include "pism/util/fem/Element.hh"
24#include "pism/util/fem/ElementIterator.hh"
25#include "pism/util/petscwrappers/SNES.hh"
26#include "pism/util/TerminationReason.hh"
27#include "pism/util/Mask.hh"
28
29namespace pism {
30
31namespace stressbalance {
32
33//! PISM's SSA solver: the finite element method implementation written by Jed and David
34/*!
35 Jed's original code is in rev 831: src/base/ssaJed/...
36 The SSAFEM duplicates most of the functionality of SSAFD, using the finite element method.
37*/
38class SSAFEM : public SSA {
39public:
40 SSAFEM(std::shared_ptr<const Grid> g);
41
42 virtual ~SSAFEM() = default;
43
44protected:
45 virtual void init_impl();
46 void cache_inputs(const Inputs &inputs);
47
48 //! Storage for SSA coefficients at element nodes.
49 //!
50 //! All fields must be "double" or structures containing "double"
51 //! for array::Array2D<T> to work correctly.
52 struct Coefficients {
53 //! ice thickness
54 double thickness;
55 //! bed elevation
56 double bed;
57 //! sea level
58 double sea_level;
59 //! basal yield stress
60 double tauc;
61 //! ice hardness
62 double hardness;
63 //! prescribed gravitational driving stress
65 };
66
69
71 double m_alpha;
72 double m_rho_g;
73
75
76 void quad_point_values(const fem::Element &E,
77 const Coefficients *x,
78 int *mask,
79 double *thickness,
80 double *tauc,
81 double *hardness) const;
82
83 static void explicit_driving_stress(const fem::Element &E,
84 const Coefficients *x,
85 Vector2d *result);
86
87 void driving_stress(const fem::Element &E, const Coefficients *x, Vector2d *result) const;
88
89 void PointwiseNuHAndBeta(double thickness,
90 double hardness,
91 int mask,
92 double tauc,
93 const Vector2d &U,
94 const Vector2d &U_x,
95 const Vector2d &U_y,
96 double *nuH, double *dnuH,
97 double *beta, double *dbeta);
98
99 void compute_local_function(Vector2d const *const * velocity,
100 Vector2d **residual);
101
102 void compute_local_jacobian(Vector2d const *const * velocity, Mat J);
103
104 virtual void solve(const Inputs &inputs);
105
106 std::shared_ptr<TerminationReason> solve_with_reason(const Inputs &inputs);
107
108 std::shared_ptr<TerminationReason> solve_nocache();
109
110 //! Adaptor for gluing SNESDAFormFunction callbacks to an SSAFEM.
111 /* The callbacks from SNES are mediated via SNESDAFormFunction, which has the
112 convention that its context argument is a pointer to a struct
113 having a DA as its first entry. The CallbackData fulfills
114 this requirement, and allows for passing the callback on to an honest
115 SSAFEM object. */
117 DM da;
119 };
120
121 // objects used internally
123
125
126 //! Storage for node types (interior, boundary, exterior).
128 //! Boundary integral (CFBC contribution to the residual).
130
134
137 // fem::P1Element m_p1_element;
138
139 // Support for direct specification of driving stress to the FEM SSA solver. This helps
140 // with certain test cases where the grid is periodic but the driving stress cannot be the
141 // gradient of a periodic function. (See commit ffb4be16.)
144private:
145 void cache_residual_cfbc(const Inputs &inputs);
146 void monitor_jacobian(Mat Jac);
147 void monitor_function(Vector2d const *const * velocity_global,
148 Vector2d const *const * residual_global);
149
150 //! SNES callbacks.
151 /*! These simply forward the call on to the SSAFEM member of the CallbackData */
152 static PetscErrorCode function_callback(DMDALocalInfo *info,
153 Vector2d const *const * velocity,
154 Vector2d **residual,
155 CallbackData *fe);
156 static PetscErrorCode jacobian_callback(DMDALocalInfo *info,
157 Vector2d const *const * velocity,
158 Mat A, Mat J, CallbackData *fe);
159};
160
161
162} // end of namespace stressbalance
163} // end of namespace pism
164
165#endif /* _SSAFEM_H_ */
This class represents a 2D vector field (such as ice velocity) at a certain grid point.
Definition Vector2d.hh:29
A storage vector combining related fields in a struct.
Definition Array2D.hh:32
Manages iterating over element indices.
The mapping from global to local degrees of freedom.
Definition Element.hh:61
Q1 element with sides parallel to X and Y axes.
Definition Element.hh:264
std::shared_ptr< TerminationReason > solve_with_reason(const Inputs &inputs)
Definition SSAFEM.cc:184
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:503
void driving_stress(const fem::Element &E, const Coefficients *x, Vector2d *result) const
Definition SSAFEM.cc:439
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:267
static PetscErrorCode function_callback(DMDALocalInfo *info, Vector2d const *const *velocity, Vector2d **residual, CallbackData *fe)
SNES callbacks.
Definition SSAFEM.cc:1211
static void explicit_driving_stress(const fem::Element &E, const Coefficients *x, Vector2d *result)
Definition SSAFEM.cc:379
CallbackData m_callback_data
Definition SSAFEM.hh:122
virtual ~SSAFEM()=default
void compute_local_jacobian(Vector2d const *const *velocity, Mat J)
Implements the callback for computing the Jacobian.
Definition SSAFEM.cc:957
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:1227
void monitor_jacobian(Mat Jac)
Definition SSAFEM.cc:1169
void monitor_function(Vector2d const *const *velocity_global, Vector2d const *const *residual_global)
Definition SSAFEM.cc:913
array::Scalar1 m_bc_mask
Definition SSAFEM.hh:67
virtual void solve(const Inputs &inputs)
Definition SSAFEM.cc:170
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:127
std::shared_ptr< TerminationReason > solve_nocache()
Definition SSAFEM.cc:194
void cache_residual_cfbc(const Inputs &inputs)
Compute and cache residual contributions from the integral over the lateral boundary.
Definition SSAFEM.cc:551
void compute_local_function(Vector2d const *const *velocity, Vector2d **residual)
Implements the callback for computing the residual.
Definition SSAFEM.cc:735
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:345
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
PISM's SSA solver: the finite element method implementation written by Jed and David.
Definition SSAFEM.hh:38
PISM's SSA solver.
Definition SSA.hh:110
const array::Vector1 & velocity() const
Get the thickness-advective 2D velocity.
static const double g
Definition exactTestP.cc:36
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