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
Blatter.hh
Go to the documentation of this file.
1/* Copyright (C) 2020, 2021, 2022, 2023, 2025 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#ifndef PISM_BLATTER_H
20#define PISM_BLATTER_H
21
22#include "pism/stressbalance/ShallowStressBalance.hh"
23#include "pism/util/petscwrappers/SNES.hh"
24#include "pism/util/petscwrappers/DM.hh"
25#include "pism/util/petscwrappers/Vec.hh"
26#include "pism/util/fem/FEM.hh"
27#include "pism/util/fem/Element.hh"
28
29namespace pism {
30
31namespace fem {
32class Element3;
33class Q1Element3Face;
34} // end of namespace fem
35
36namespace stressbalance {
37
39public:
40 Blatter(std::shared_ptr<const Grid> grid, int Mz, int coarsening_factor);
41 virtual ~Blatter() = default;
42
43 void update(const Inputs &inputs, bool /*full_update*/);
44
45 std::shared_ptr<array::Array3D> velocity_u_sigma() const;
46 std::shared_ptr<array::Array3D> velocity_v_sigma() const;
47
48 /*!
49 * 2D input parameters
50 */
51 struct Parameters {
52 // elevation (z coordinate) of the bottom domain boundary
53 double bed;
54 // thickness of the domain
55 double thickness;
56 // NodeType stored as double
57 double node_type;
58 // basal yield stress
59 double tauc;
60 // sea level elevation (used to determine if a location is grounded)
61 double sea_level;
62 // floatation function (positive where floating, zero or negative where grounded)
63 double floatation;
64 // FIXME: Add Dirichlet BC at a map plane location.
65 };
66
67protected:
68 // u and v components of ice velocity on the fine sigma grid
69 std::shared_ptr<array::Array3D> m_u_sigma, m_v_sigma;
70
71 // 3D dof=2 DM used by m_snes
73 // storage for the solution
75 // storage for the old solution during parameter continuation
77 // solver
79
81
82 // Scaling of quadrature weights (note: this does not seem to matter).
83 double m_scaling;
84
85 // Ice density times g
87
88 // Water density times g
90
91 // The flow law Glen exponent
93
94 // Use the eta-transformation to compute surface gradient
96
97 // Viscosity regularization constant
99
100 // Enhancement factor for ice viscosity.
102
103 // True if the Eisenstat-Walker method of adjusting linear solver tolerances is enabled.
105
106 static const int m_Nq = 100;
107 static const int m_n_work = 9;
108
110
112
115
116 void init_impl();
117
118 void define_model_state_impl(const File &output) const;
119
120 void write_model_state_impl(const File &output) const;
121
122 static bool exterior_element(const int *node_type);
123
124 static bool grounding_line(const double *F);
125
126 static bool partially_submerged_face(int face, const double *z, const double *sea_level);
127
128 void compute_node_type(double min_thickness);
129
130 virtual void nodal_parameter_values(const fem::Q1Element3 &element,
131 Parameters **P,
132 int i,
133 int j,
134 int *node_type,
135 double *bottom,
136 double *thickness,
137 double *surface,
138 double *sea_level) const;
139
140 virtual bool marine_boundary(int face,
141 const int *node_type,
142 const double *ice_bottom,
143 const double *sea_level);
144
145 virtual bool dirichlet_node(const DMDALocalInfo &info, const fem::Element3::GlobalIndex& I);
146
147 virtual Vector2d u_bc(double x, double y, double z) const;
148
149 void compute_jacobian(DMDALocalInfo *info, const Vector2d ***x, Mat A, Mat J);
150
151 void jacobian_dirichlet(const DMDALocalInfo &info, Parameters **P, Mat J);
152
153 virtual void jacobian_f(const fem::Q1Element3 &element,
154 const Vector2d *u_nodal,
155 const double *B_nodal,
156 double K[2 * fem::q13d::n_chi][2 * fem::q13d::n_chi]);
157
158 virtual void jacobian_basal(const fem::Q1Element3Face &face,
159 const double *tauc_nodal,
160 const double *f_nodal,
161 const Vector2d *u_nodal,
162 double K[2 * fem::q13d::n_chi][2 * fem::q13d::n_chi]);
163
164 void compute_residual(DMDALocalInfo *info, const Vector2d ***X, Vector2d ***R);
165
166 void residual_dirichlet(const DMDALocalInfo &info,
167 Parameters **P,
168 const Vector2d ***x,
169 Vector2d ***R);
170
171 virtual void residual_f(const fem::Q1Element3 &element,
172 const Vector2d *u_nodal,
173 const double *B_nodal,
174 Vector2d *residual);
175
176 virtual void residual_source_term(const fem::Q1Element3 &element,
177 const double *surface,
178 const double *bed,
179 Vector2d *residual);
180
181 virtual void residual_basal(const fem::Q1Element3 &element,
182 const fem::Q1Element3Face &face,
183 const double *tauc_nodal,
184 const double *f_nodal,
185 const Vector2d *u_nodal,
186 Vector2d *residual);
187
188 virtual void residual_surface(const fem::Q1Element3 &element,
189 const fem::Q1Element3Face &face,
190 Vector2d *residual);
191
192 virtual void residual_lateral(const fem::Q1Element3 &element,
193 const fem::Q1Element3Face &face,
194 const double *surface_nodal,
195 const double *z_nodal,
196 const double *sl_nodal,
197 Vector2d *residual);
198
199 static PetscErrorCode jacobian_callback(DMDALocalInfo *info,
200 const Vector2d ***x,
201 Mat A, Mat J,
202 Blatter *solver);
203
204 static PetscErrorCode function_callback(DMDALocalInfo *info,
205 const Vector2d ***x, Vector2d ***f,
206 Blatter *solver);
207
208 virtual void init_2d_parameters(const Inputs &inputs);
209
210 void init_ice_hardness(const Inputs &inputs, const petsc::DM &da);
211
212 // Guts of the constructor. This method wraps PETSc calls to simplify error checking.
213 PetscErrorCode setup(DM pism_da, grid::Periodicity p, int Mz, int coarsening_factor,
214 const std::string &prefix);
215
216 void set_initial_guess(const array::Array3D &u_sigma, const array::Array3D &v_sigma);
217
218 void copy_solution();
219
221
222 void get_basal_velocity(array::Vector &result);
223
224 void report_mesh_info();
225
228 SNESConvergedReason snes_reason;
229 KSPConvergedReason ksp_reason;
230 };
231
234};
235
236} // end of namespace stressbalance
237} // end of namespace pism
238
239#endif /* PISM_BLATTER_H */
std::shared_ptr< const Grid > grid() const
Definition Component.cc:105
High-level PISM I/O class.
Definition File.hh:55
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
A virtual class collecting methods common to ice and bedrock 3D fields.
Definition Array3D.hh:33
3D Q1 element
Definition Element.hh:358
void jacobian_dirichlet(const DMDALocalInfo &info, Parameters **P, Mat J)
Definition jacobian.cc:169
void define_model_state_impl(const File &output) const
The default (empty implementation).
Definition Blatter.cc:692
void get_basal_velocity(array::Vector &result)
Definition Blatter.cc:1127
virtual void residual_f(const fem::Q1Element3 &element, const Vector2d *u_nodal, const double *B_nodal, Vector2d *residual)
Definition residual.cc:37
virtual Vector2d u_bc(double x, double y, double z) const
Definition Blatter.cc:133
static const int m_n_work
Definition Blatter.hh:107
virtual void residual_source_term(const fem::Q1Element3 &element, const double *surface, const double *bed, Vector2d *residual)
Definition residual.cc:95
void compute_node_type(double min_thickness)
Definition Blatter.cc:57
virtual void residual_basal(const fem::Q1Element3 &element, const fem::Q1Element3Face &face, const double *tauc_nodal, const double *f_nodal, const Vector2d *u_nodal, Vector2d *residual)
Definition residual.cc:164
void write_model_state_impl(const File &output) const
The default (empty implementation).
Definition Blatter.cc:697
void compute_residual(DMDALocalInfo *info, const Vector2d ***X, Vector2d ***R)
Definition residual.cc:312
virtual bool marine_boundary(int face, const int *node_type, const double *ice_bottom, const double *sea_level)
Definition Blatter.cc:226
fem::Q1Element3Face m_face4
Definition Blatter.hh:113
static bool partially_submerged_face(int face, const double *z, const double *sea_level)
Definition Blatter.cc:193
virtual void residual_lateral(const fem::Q1Element3 &element, const fem::Q1Element3Face &face, const double *surface_nodal, const double *z_nodal, const double *sl_nodal, Vector2d *residual)
Definition residual.cc:218
virtual void init_2d_parameters(const Inputs &inputs)
Definition Blatter.cc:533
std::shared_ptr< array::Array3D > m_u_sigma
Definition Blatter.hh:69
std::shared_ptr< array::Array3D > velocity_v_sigma() const
Definition Blatter.cc:1206
SolutionInfo parameter_continuation()
Definition Blatter.cc:816
void update(const Inputs &inputs, bool)
Definition Blatter.cc:983
PetscErrorCode setup(DM pism_da, grid::Periodicity p, int Mz, int coarsening_factor, const std::string &prefix)
Definition Blatter.cc:383
double m_work[m_n_work][m_Nq]
Definition Blatter.hh:109
void init_ice_hardness(const Inputs &inputs, const petsc::DM &da)
Definition Blatter.cc:580
static PetscErrorCode jacobian_callback(DMDALocalInfo *info, const Vector2d ***x, Mat A, Mat J, Blatter *solver)
Definition jacobian.cc:355
void residual_dirichlet(const DMDALocalInfo &info, Parameters **P, const Vector2d ***x, Vector2d ***R)
Definition residual.cc:263
static PetscErrorCode function_callback(DMDALocalInfo *info, const Vector2d ***x, Vector2d ***f, Blatter *solver)
Definition residual.cc:458
static const int m_Nq
Definition Blatter.hh:106
fem::Q1Element3Face m_face100
Definition Blatter.hh:114
void compute_averaged_velocity(array::Vector &result)
Definition Blatter.cc:1167
std::shared_ptr< array::Array3D > velocity_u_sigma() const
Definition Blatter.cc:1202
static bool exterior_element(const int *node_type)
Definition Blatter.cc:147
virtual void jacobian_f(const fem::Q1Element3 &element, const Vector2d *u_nodal, const double *B_nodal, double K[2 *fem::q13d::n_chi][2 *fem::q13d::n_chi])
Definition jacobian.cc:36
static bool grounding_line(const double *F)
Definition Blatter.cc:167
virtual void jacobian_basal(const fem::Q1Element3Face &face, const double *tauc_nodal, const double *f_nodal, const Vector2d *u_nodal, double K[2 *fem::q13d::n_chi][2 *fem::q13d::n_chi])
Definition jacobian.cc:123
void compute_jacobian(DMDALocalInfo *info, const Vector2d ***x, Mat A, Mat J)
Definition jacobian.cc:200
array::Array2D< Parameters > m_parameters
Definition Blatter.hh:80
virtual bool dirichlet_node(const DMDALocalInfo &info, const fem::Element3::GlobalIndex &I)
Definition Blatter.cc:125
virtual void residual_surface(const fem::Q1Element3 &element, const fem::Q1Element3Face &face, Vector2d *residual)
Definition residual.cc:201
std::shared_ptr< array::Array3D > m_v_sigma
Definition Blatter.hh:69
virtual ~Blatter()=default
Vector2d m_work2[m_n_work][m_Nq]
Definition Blatter.hh:111
virtual void nodal_parameter_values(const fem::Q1Element3 &element, Parameters **P, int i, int j, int *node_type, double *bottom, double *thickness, double *surface, double *sea_level) const
Definition Blatter.cc:636
void set_initial_guess(const array::Array3D &u_sigma, const array::Array3D &v_sigma)
Definition Blatter.cc:1143
Shallow stress balance (such as the SSA).
const int n_chi
Number of shape functions on a Q1 element.
Definition FEM.hh:213
Periodicity
Definition Grid.hh:54
static double F(double SL, double B, double H, double alpha)