Blatter stress balance solver: technical details¶
Notation¶
2D vector |
|
outward-pointing normal vector at domain boundaries |
|
ice viscosity (see (102)) |
|
ice density |
|
water (usually sea water) density |
|
gravitational acceleration |
|
sea level elevation |
|
ice thickness |
|
ice surface elevation |
|
ice hardness |
|
ice softness |
|
regularization parameter |
|
basal resistance coefficient |
|
Glen exponent |
Introduction¶
This implementation is based on the PETSc example snes/tutorials/ex48.c
(see
[56]) and is inspired by [175],
[73], [75], [72], and [176].
Define
Using this notation, the Blatter-Pattyn (BP) stress balance equations are
Here “
where
It is also assumed that
ice is incompressible and
, and and .
The BP approximation of the second invariant
Weak form¶
Note
Recall the product rule
and the divergence theorem:
We omit discussions of function spaces; see [175], [176] and other references mentioned above for details.
To obtain the weak form, we multiply both equations in (101) by a scalar test
function
Similarly, multiplying the second equation by
We combine these and say that the weak form of (101) is
where
The first term corresponds to Neumann and Robin boundary conditions; it
vanishes for “natural” BC
Boundary conditions¶
The domain boundary consists of three disjoint parts:
The interface between ice and the underlying bed or (if the ice is floating) water.
At this interface PISM uses Robin BC implementing basal sliding.
The interface between ice and the air[1] above it.
The integral over this part of the boundary is zero because we assume that natural (“no stress”) boundary conditions apply. (We ignore atmospheric pressure.)
Vertical “cliffs” at grounded margins are interpreted as approximations of the very steep, but not vertical, upper surface. Following this interpretation we use natural boundary conditions at grounded lateral margins.
The vertical interface between ice and air (above sea level) or water (below sea level) at marine ice margins.
At this interface PISM uses Neumann BC corresponding to the difference between the cryostatic pressure of the ice on one side and the hydrostatic pressure of water on the other side of the interface.
In addition to this we support Dirichlet boundary conditions for verification and to de-couple unknowns at ice-free locations.
Note
Our implementation supports Dirichlet boundary conditions, but this feature is not exposed to the rest of PISM.
Unlike [56] we do not support “no-slip” BC at the base. This allows us to avoid Jacobian scaling tricks they needed to achieve good multigrid performance.
For each node that belongs to the Dirichlet boundary we assemble “trivial” equations
where
The implementation avoids adding contributions from adjacent elements to residual and Jacobian entries corresponding to Dirichlet locations. Prescribed velocity values are substituted into equations that depend on them (see section 3.2 of [56] for details).
Basal boundary¶
The boundary condition corresponding to basal sliding is
In the weak form (106) this corresponds to replacing the first term:
Here
Where ice is grounded
The grounding line (if present) divides bottom faces of some elements into grounded parts
that experience drag and floating parts that do not. This implementation uses a low-order
quadrature with many equally-spaced points (a Newton-Cotes quadrature) to integrate over
the part of the basal boundary containing the grounding line. Here
Note
It may be a good idea to implement scaling of the
coefficient according to the fraction of an element face that is grounded, similar to SEP1 (equation 7 in [178]). An approximation of the grounded fraction for an element column can be computed using existing code in PISM.In general the bottom face of an element is not planar and not parallel to the plane
. This means that integrals over the basal boundary should be evaluated using parameterizations of element faces (i.e. as surface integrals) and not using 2D FEM machinery.
Marine margins¶
We assume that marine ice margins consist of vertical cliffs, i.e. the outward-pointing
normal vector has the form
The Neumann boundary condition at marine margins is
In other words, this boundary condition corresponds to the difference between the cryostatic pressure of the ice on one side of the interface and the hydrostatic pressure of water on the other. The atmospheric pressure is ignored. Equation (109) is a generalization of equation 18 in [175].
Just like in the implementation of the basal boundary condition near grounding lines, we use a low order quadrature with many equally-spaced points to approximate integrals over element faces intersected by the sea level. This should improve the quality of the approximation of this boundary condition (note that the right hand side of (109) is continuous but not continuously differentiable).
Note
We need to evaluate the importance of the quadrature choice described above. Does using depth-dependent BC matter? (We could simplify the code and use depth-averaged BC it if it does not.) Should we use lots of quadrature points?
Note
CISM 2.1 [73] uses a depth-averaged lateral boundary condition without justification.
The lateral BC described in [175] is equivalent to the one described here, but the implementation in Albany/FELIX (and therefore in MALI) matches the one in CISM.
Implementations in CISM and MALI do not use this boundary condition at grounded margins because doing so appears to produce over-estimates of the ice speed near grounded ice margins. Our experiments show the same behavior.
Solver implementation¶
Discretization¶
To create a non-linear algebraic system of equations approximating (106),
we create a hexahedral mesh on the domain
Let
Then the problem is
Find
, ( ) so that
holds for all
As in section 3 of [56], we write the discretization of
(106) as an algebraic system
where
(compare to (106)) and
Residual evaluation¶
The residual evaluation is performed in the usual manner, by iterating over all the elements containing ice (see Model domain and mesh structure) and adding element contributions to the “global” residual vector.
The residual itself can be broken up into the following parts:
The basal boundary term implementing the basal drag (part of
)The “main” part (
)The source term corresponding to the driving stress (
)The top boundary part (zero in actual simulations because we use the natural BC at the top boundary; can be non-zero in verification tests, part of
)The marine boundary part implementing stress (Neumann) BC at the calving front (part of
).Residual at Dirichlet nodes.
This decomposition makes it possible to use source terms dictated by the choice of a manufactured solution while keeping (and testing) the rest of the code.
Note
We integrate over the whole domain (
Main residual contribution¶
Here
Driving stress contribution¶
Basal contribution¶
Marine boundary contribution¶
where
Residual at Dirichlet BC locations¶
where
Jacobian evaluation¶
We use an analytical (as opposed to approximated using finite differences or automatic differentiation) Jacobian computed using formulas listed below.
Here we focus on the derivation of the Jacobian contribution corresponding to the “main”
part of the residual (
This Jacobian contribution has four parts:
with
Let
and (using the product rule) we get
The derivatives of
Taking derivatives of
and
The derivatives of
To compute
and so on.
Basal contribution to the Jacobian¶
Here
Note that
Recall (see (115)) that the basal sliding contribution to the residual is equal to
and so the basal contribution to the Jacobian consists of partial derivatives of these
with respect to
For example,
Jacobian at Dirichlet locations¶
The Jacobian at Dirichlet locations is set to
More specifically (due to the interleaved ordering of the unknowns:
Note
It may be interesting to see if varying the scaling of Jacobian entries at Dirichlet
nodes affects the condition number of the Jacobian matrix. See the variable scaling
in the code and set
-bp_ksp_view_singularvalues
to see if it has an effect.
Iceberg elimination¶
As described in [72], isolated patches of ice with low basal resistance and patches connected only via a single node (a “hinge”) are problematic because the system (101) determines ice velocity up to rigid rotations and translations.
This is not a new issue: both FD and FEM solvers of the SSA stress balance require some form of iceberg elimination. We use a connected component labeling algorithm to identify patches of floating ice. In the FEM context this requires inspecting elements: two elements are considered connected if they share a boundary.
Note
We could improve this mechanism by implementing a version of the method described in [72]: instead of removing floating ice not connected to grounded ice (PISM’s approach) they
Identify all the connected patches of ice.
Remove patches of ice which have no mesh nodes with
Preconditioning¶
Back in 2013 Brown et al [56] showed that multigrid can act as an effective preconditioner for BP stress balance solvers. However, all test cases considered in that work use periodic boundary conditions and therefore avoid all the issues related to the moving ice margin present in complete ice sheet models. Our tests suggests that a naive implementation assembling trivial equations at “ice free” nodes combined with standard geometric multigrid (coarsening in all 3 directions) is not likely to succeed and we need a different approach.
So, we use semi-coarsening in the vertical direction even though Brown et al state that “semi-coarsening is unattractive”. One of the arguments against semi-coarsening is the larger number of multigrid levels needed: semi-coarsening gives a smaller reduction in the number of unknowns from one level to the next (factor of 2 instead of 8 in the full multigrid approach). Our tests show that “aggressive” semi-coarsening (i.e. using coarsening factors larger than 2 and as high as 8) appears to be effective, allowing one to achieve similar reductions in the number of unknowns from one level to the next in a multigrid hierarchy.
The second argument against semi-coarsening is deeper: spatial variations in the sliding
parameter
Inspired by [72], we use geometric multigrid to build a mesh hierarchy
with the coarsest level containing a small number (2 or 3) of vertical levels combined
with algebraic multigrid as a preconditioner for the solver on the coarsest level.
(Semi-coarsening in the vertical direction cannot reduce the number of unknowns in the
In addition to this, we follow [56] in ordering unknowns so that
columns are contiguous (and
Vertical grid sizes compatible with coarsening¶
Ideally, the coarsest mesh in the hierarchy should have 2 nodes in the
for some positive integer
nodes in the
This means that for a given stress_balance
.blatter
.coarsening_factor
and number
of multigrid levels -bp_pc_mg_levels k
the value of
Controlling using PETSc options¶
The PETSc SNES
object solving the BP system uses the bp_
prefix for all
command-line options.
To choose a preconditioner, set
-bp_pc_type XXX
where XXX
is the name of a preconditioner.
Run
pism -stress_balance blatter [other options] -help | grep "-bp"
to get the list of all PETSc options controlling this solver.
To use a geometric multigrid preconditioner with
-bp_pc_type mg -bp_pc_mg_levels N
An “aggressive” (i.e. greater than 2) coarsening factor may work well. Use
stress_balance
.blatter
.coarsening_factor
to set it.
See Vertical grid sizes compatible with coarsening for the discussion of the relationship between the number of vertical levels, number of multigrid levels, and the coarsening factor.
Set
-bp_mg_coarse_pc_type gamg
to use PETSc’s GAMG on the coarsest multigrid level.
Note
It would be interesting to compare different preconditioning options on the coarsest MG level (GAMG, Hypre BoomerAMG, …).
Additional code needed to support geometric multigrid¶
To support the multigrid preconditioner we need to be able re-discretize the system on the mesh provided by PETSc in our to residual and Jacobian evaluators. In general, this requires
re-computing grid-related constants (
, etc) using the grid (avoid using hard-wired constants, e.g. computed using the fine grid), andadditional code to restrict gridded inputs from the fine grid mesh to coarser meshes.
This solver does not support coarsening in horizontal directions, so gridded
two-dimensional inputs can be used on all multigrid levels. The grid spacing (
To transfer the one three-dimensional gridded input field (ice hardness), we create interpolation matrices mapping from a coarse level to the next (finer) level in the hierarchy. The transpose of this matrix is used as a restriction operator.
Parameter continuation as a recovery mechanism¶
As in [175], we can start with a large
Note
Not implemented yet.
Model domain and mesh structure¶
The domain is
where
Coordinates of the mesh nodes have the form
Each element’s projection onto the plane
Supporting evolving ice extent¶
The ice extent changes as the model runs and the solver implementation has to allow for these changes.
To simplify the logic used to identify the interior of the ice volume and its lateral
boundaries we compute the “type” of each node in the mesh. Given a threshold stress_balance
.ice_free_thickness_standard
), we say that
an element contains ice if ice thickness at all its nodes equals to or exceeds
,a node is interior (within the ice) if all the elements it belongs to contain ice,
a node is exterior (outside the ice volume) if no element it belongs to contains ice,
a node that is neither interior nor exterior is a boundary node.
Only elements containing ice are included in the residual and Jacobian evaluation.
We prescribe zero
In addition to this, we need to identify vertical faces of elements at lateral boundaries.
An element face is a part of the lateral boundary if all four of its nodes are boundary nodes.
Solver inputs and outputs¶
The BP solver uses the following inputs:
basal yield stress (
),ice thickness,
bed elevation,
sea level elevation,
ice enthalpy (used to compute ice hardness)
And provides the following outputs:
and components of ice velocity at the nodes of the mesh (saved to output files to be used as an initial guess later)vertically-averaged
and (used in the mass continuity step, i.e. to update ice geometry)basal frictional heating (used to couple stress balance to energy conservation)
Steps performed by the solver¶
Compute ice bottom elevation.
Compute floatation function
( if ice is grounded, if floating).Compute node type.
Compute ice hardness at the nodes of the mesh.
Call PETSc’s
SNESSolve()
.Extract basal velocity and compute basal frictional heating.
Compute vertically-averaged ice velocity.
Integration with the rest of PISM¶
Conservation of energy¶
Coupling to PISM’s energy balance models requires
all 3 components of ice velocity on PISM’s grid, and
strain heating.
These are computed by using piecewise-linear interpolation in the vertical direction to
put
Testing and verification¶
Verification test XY¶
- Exact solution:
- Domain:
,- Boundary conditions:
Dirichlet BC corresponding to the exact solution on all lateral boundaries (
, , , ). Natural BC at the top and bottom boundaries.- Comments:
This test uses a constant ice hardness,
, and has no variation in the direction. It is similar to the MMS verification test in [175], section 4.1 (we use Dirichlet boundary conditions along the whole lateral boundary instead of Robin conditions derived from the chosen exact solution).The compensatory term is computed using SymPy; please see the code in
src/stressbalance/blatter/verification
.
Verification test XZ¶
- Exact solution:
- Domain:
- Boundary conditions:
Dirichlet BC corresponding to the exact solution at
and .Basal (
) BC is a combination (i.e. sum) of the BC in Basal boundary and a compensatory term derived using the exact solution.Top surface (
) BC is derived using the exact solution.
- Comments:
This test
uses a constant ice hardness,
Glen exponent
,has no variation (and is periodic) in the
direction,uses a constant basal resistance coefficient
.
It is similar to the
MMS verification test [175], section 4.2 (again, we use Dirichlet BC at lateral boundaries instead of Robin conditions stated in the paper).See [175] and the code in
src/stressbalance/blatter/verification
for details.
Verification test XZ-CFBC¶
This setup tests the “calving front boundary condition” (see Marine margins).
- Exact solution:
- Domain:
,- Boundary conditions:
Dirichlet BC at
.Uses the BC described in Marine margins at
.
- Comments:
This test uses the Glen exponent of
(constant viscosity) and has no variation in the direction.The sea level is set to
, overriding the floatation criterion to ensure that all the ice is submerged.
Verification test XZ-VV (van der Veen profile)¶
This setup tests the implementation of the basal sliding boundary condition (see Basal boundary) using the van der Veen shelf profile [180]:
where
The surface elevation
for some positive constant
Note
Functions
if the Glen exponent
We use a Dirichlet BC at the left boundary:
and a stress BC at the right boundary:
The stress BC at the top surface is
The boundary condition at the bottom surface has the form
Known issues and future work¶
Eliminate “wiggles” near areas with steep surface slopes.
Implement drag at lateral boundaries in fjords and alpine valleys.
Implement parameter continuation.
Couple to melange back pressure parameterizations by replacing
in Marine margins.
Footnotes
Previous | Up | Next |