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"
32#include "pism/util/node_types.hh"
33#include "pism/util/pism_utilities.hh"
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"
41namespace stressbalance {
50 m_bc_mask(grid,
"bc_mask"),
51 m_bc_values(grid,
"_bc"),
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()) {
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");
82 ierr = DMDASNESSetFunctionLocal(*dm, INSERT_VALUES,
83#
if PETSC_VERSION_LT(3,21,0)
89 PISM_CHK(ierr,
"DMDASNESSetFunctionLocal");
91 ierr = DMDASNESSetJacobianLocal(*dm,
92#
if PETSC_VERSION_LT(3,21,0)
98 PISM_CHK(ierr,
"DMDASNESSetJacobianLocal");
100 ierr = DMSetMatType(*dm,
"baij");
104 PISM_CHK(ierr,
"DMSetApplicationContext");
106 ierr = SNESSetDM(
m_snes, *dm);
110 int snes_max_it = 200;
111 ierr = SNESSetTolerances(
m_snes, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, snes_max_it,
113 PISM_CHK(ierr,
"SNESSetTolerances");
115 ierr = SNESSetFromOptions(
m_snes);
116 PISM_CHK(ierr,
"SNESSetFromOptions");
119 "node types: interior, boundary, exterior");
124 "residual contribution from lateral boundaries");
133 if (
m_grid->variables().is_available(
"ssa_driving_stress_x") and
134 m_grid->variables().is_available(
"ssa_driving_stress_y")) {
139 m_log->message(2,
" [using the SNES-based finite element method implementation]\n");
145 "Enforce Dirichlet conditions with this additional scaling",
174 if (reason->failed()) {
176 "SSAFEM solve failed to converge (SNES reason %s)",
177 reason->description().c_str());
180 if (
m_log->get_threshold() > 2) {
181 m_stdout_ssa +=
"SSAFEM converged (SNES reason " + reason->description() +
")";
203 ierr = PetscViewerASCIIOpen(
m_grid->com, filename->c_str(), viewer.
rawptr());
204 PISM_CHK(ierr,
"PetscViewerASCIIOpen");
206 ierr = PetscViewerASCIIPrintf(viewer,
"SNES before SSASolve_FE\n");
207 PISM_CHK(ierr,
"PetscViewerASCIIPrintf");
209 ierr = SNESView(
m_snes, viewer);
212 ierr = PetscViewerASCIIPrintf(viewer,
"solution vector before SSASolve_FE\n");
213 PISM_CHK(ierr,
"PetscViewerASCIIPrintf");
220 if (
m_log->get_threshold() >= 2) {
229 SNESConvergedReason snes_reason;
230 ierr = SNESGetConvergedReason(
m_snes, &snes_reason);
PISM_CHK(ierr,
"SNESGetConvergedReason");
233 if (not reason->failed()) {
238 bool view_solution =
options::Bool(
"-ssa_view_solution",
"view solution of the SSA system");
241 ierr = PetscViewerASCIIOpen(
m_grid->com, filename->c_str(), viewer.
rawptr());
242 PISM_CHK(ierr,
"PetscViewerASCIIOpen");
244 ierr = PetscViewerASCIIPrintf(viewer,
"solution vector after SSASolve\n");
245 PISM_CHK(ierr,
"PetscViewerASCIIPrintf");
279 const std::vector<double> &z =
m_grid->z();
289 if (use_explicit_driving_stress) {
295 for (
auto p :
m_grid->points()) {
296 const int i = p.i(), j = p.j();
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);
311 m_grid->kBelowHeight(thickness),
331 const bool use_cfbc =
m_config->get_flag(
"stress_balance.calving_front_stress_bc");
335 m_config->get_number(
"stress_balance.ice_free_thickness_standard"),
351 double *hardness)
const {
353 unsigned int n = E.
n_pts();
355 for (
unsigned int q = 0; q <
n; q++) {
364 for (
int k = 0;
k < E.
n_chi();
k++) {
374 mask[q] =
m_gc.
mask(sea_level, bed, thickness[q]);
383 const unsigned int n = E.
n_pts();
385 for (
unsigned int q = 0; q <
n; q++) {
388 for (
int k = 0;
k < E.
n_chi();
k++) {
443 const unsigned int n = E.
n_pts();
445 for (
unsigned int q = 0; q <
n; q++) {
457 for (
int k = 0;
k < E.
n_chi();
k++) {
471 const int M =
m_gc.
mask(sea_level, b, 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;
479 result[q].
u = - pressure * h_x;
480 result[q].
v = - pressure * h_y;
511 double *nuH,
double *dnuH,
512 double *beta,
double *dbeta) {
514 if (thickness < strength_extension->get_min_thickness()) {
516 if (dnuH !=
nullptr) {
525 if (dnuH !=
nullptr) {
539 if (dbeta !=
nullptr) {
559 const unsigned int Nq = 2;
560 double chi_b[Nq][Nq];
565 for (
int i = 0; i < 2; ++i) {
566 for (
int j = 0; j < 2; ++j) {
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)};
584 use_cfbc =
m_config->get_flag(
"stress_balance.calving_front_stress_bc");
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");
614 for (
int j = ys; j < ys + ym; j++) {
615 for (
int i = xs; i < xs + xm; i++) {
633 E = &p1_element[type];
643 std::vector<Vector2d> I(Nk);
655 double psi[2] = {0.0, 0.0};
657 const unsigned int n_sides = E->
n_sides();
659 for (
unsigned int s = 0; s < n_sides; ++s) {
667 n1 = (n0 + 1) % n_sides;
675 for (
unsigned int q = 0; q < Nq; ++q) {
681 psi[0] = chi_b[0][q];
682 psi[1] = chi_b[1][q];
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];
693 P_ice = 0.5 * ice_density * standard_gravity * H,
695 ice_density, ocean_density,
699 double dP = H * (P_ice - P_water);
711 I[n0] += W * (- psi[0] * dP) * E->
normal(s);
712 I[n1] += W * (- psi[1] * dP) * E->
normal(s);
741 const bool use_cfbc =
m_config->get_flag(
"stress_balance.calving_front_stress_bc");
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)};
758 for (
auto p :
m_grid->points()) {
759 const int i = p.i(), j = p.j();
768 Vector2d U[Nq_max], U_x[Nq_max], U_y[Nq_max];
779 for (
int j = ys; j < ys + ym; j++) {
780 for (
int i = xs; i < xs + xm; i++) {
800 E = &p1_element[type];
811 const unsigned int Nq = E->
n_pts();
817 double thickness[Nq_max];
819 double hardness[Nq_max];
828 if (use_explicit_driving_stress) {
841 if (dirichlet_data) {
844 dirichlet_data.
enforce(*E, velocity_nodal);
856 for (
unsigned int k = 0;
k < Nk;
k++) {
861 for (
unsigned int q = 0; q < Nq; q++) {
865 double eta = 0.0, beta = 0.0;
867 U[q], U_x[q], U_y[q],
868 &eta, NULL, &beta, NULL);
871 const Vector2d tau_b = U[q] * (- beta);
876 u_y_plus_v_x = U_y[q].
u + U_x[q].
v;
879 for (
int k = 0;
k < E->
n_chi();
k++) {
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));
899 if (dirichlet_data) {
900 dirichlet_data.
fix_residual(velocity_global, residual_global);
915 Vector2d const *
const *
const residual_global) {
917 bool monitorFunction =
options::Bool(
"-ssa_monitor_function",
"monitor the SSA residual");
918 if (not monitorFunction) {
922 ierr = PetscPrintf(
m_grid->com,
923 "SSA Solution and Function values (pointwise residuals)\n");
928 for (
auto p :
m_grid->points()) {
929 const int i = p.i(), j = p.j();
931 ierr = PetscSynchronizedPrintf(
m_grid->com,
932 "[%2d, %2d] u=(%12.10e, %12.10e) f=(%12.4e, %12.4e)\n",
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");
943 ierr = PetscSynchronizedFlush(
m_grid->com, NULL);
944 PISM_CHK(ierr,
"PetscSynchronizedFlush");
960 const bool use_cfbc =
m_config->get_flag(
"stress_balance.calving_front_stress_bc");
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)};
973 PetscErrorCode ierr = MatZeroEntries(Jac);
982 Vector2d U[Nq_max], U_x[Nq_max], U_y[Nq_max];
993 for (
int j = ys; j < ys + ym; j++) {
994 for (
int i = xs; i < xs + xm; i++) {
1014 E = &p1_element[type];
1030 double thickness[Nq_max];
1031 double tauc[Nq_max];
1032 double hardness[Nq_max];
1039 mask, thickness, tauc, hardness);
1049 if (dirichlet_data) {
1050 dirichlet_data.
enforce(*E, velocity_nodal);
1054 E->
evaluate(velocity_nodal, U, U_x, U_y);
1060 double K[2*Nk][2*Nk];
1062 ierr = PetscMemzero(K,
sizeof(K));
1065 for (
unsigned int q = 0; q < Nq; q++) {
1073 u_y_plus_v_x = U_y[q].
u + U_x[q].
v;
1075 double eta = 0.0, deta = 0.0, beta = 0.0, dbeta = 0.0;
1077 U[q], U_x[q], U_y[q],
1078 &eta, &deta, &beta, &dbeta);
1080 for (
unsigned int l = 0; l < n_chi; l++) {
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;
1092 eta_u = deta * gamma_u,
1093 eta_v = deta * gamma_v;
1097 taub_xu = -dbeta * u * u *
phi.val - beta *
phi.val,
1098 taub_xv = -dbeta * u * v *
phi.val,
1099 taub_yu = -dbeta * v * u *
phi.val,
1100 taub_yv = -dbeta * v * v *
phi.val - beta *
phi.val;
1102 for (
unsigned int k = 0;
k < n_chi;
k++) {
1109 ierr = PetscPrintf(PETSC_COMM_SELF,
"eta=0 i %d j %d q %d k %d\n", i, j, q,
k);
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);
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);
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);
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);
1143 if (dirichlet_data) {
1155 ierr = MatAssemblyBegin(Jac, MAT_FINAL_ASSEMBLY);
1156 PISM_CHK(ierr,
"MatAssemblyBegin");
1158 ierr = MatAssemblyEnd(Jac, MAT_FINAL_ASSEMBLY);
1161 ierr = MatSetOption(Jac, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_TRUE);
1164 ierr = MatSetOption(Jac, MAT_SYMMETRIC, PETSC_TRUE);
1171 PetscErrorCode ierr;
1172 bool mon_jac =
options::Bool(
"-ssa_monitor_jacobian",
"monitor the SSA Jacobian");
1181 ierr = SNESGetIterationNumber(
m_snes, &iter);
1182 PISM_CHK(ierr,
"SNESGetIterationNumber");
1184 auto file_name =
pism::printf(
"PISM_SSAFEM_J%d.m", (
int)iter);
1187 "writing Matlab-readable file for SSAFEM system A xsoln = rhs to file `%s' ...\n",
1192 ierr = PetscViewerSetType(viewer, PETSCVIEWERASCII);
1193 PISM_CHK(ierr,
"PetscViewerSetType");
1195 ierr = PetscViewerPushFormat(viewer, PETSC_VIEWER_ASCII_MATLAB);
1196 PISM_CHK(ierr,
"PetscViewerPushFormat");
1198 ierr = PetscViewerFileSetName(viewer, file_name.c_str());
1199 PISM_CHK(ierr,
"PetscViewerFileSetName");
1201 ierr = PetscObjectSetName((PetscObject) Jac,
"A");
1202 PISM_CHK(ierr,
"PetscObjectSetName");
1204 ierr = MatView(Jac, viewer);
1207 ierr = PetscViewerPopFormat(viewer);
1208 PISM_CHK(ierr,
"PetscViewerPopFormat");
1213 Vector2d const *
const *
const velocity,
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");
1229 Vector2d const *
const *
const velocity,
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");
const units::System::Ptr m_sys
unit system used by this component
std::shared_ptr< const Config > m_config
configuration database used by this component
const std::shared_ptr< const Grid > m_grid
grid used by this component
std::shared_ptr< const Logger > m_log
logger (for easy access)
int mask(double sea_level, double bed, double thickness) const
array::Scalar1 sea_level_elevation
array::Scalar2 ice_thickness
array::Scalar2 bed_elevation
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
This class represents a 2D vector field (such as ice velocity) at a certain grid point.
Makes sure that we call begin_access() and end_access() for all accessed array::Arrays.
void copy_from(const Array2D< T > &source)
void add(double alpha, const Array2D< T > &x)
double * get_column(int i, int j)
void set_interpolation_type(InterpolationType type)
void set(double c)
Result: v[j] <- c for all j.
std::shared_ptr< petsc::DM > dm() const
VariableMetadata & metadata(unsigned int N=0)
Returns a reference to the VariableMetadata object containing metadata for the compoment N.
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.
void reset(int i, int j)
Initialize the Element to element (i, j) for the purposes of inserting into global residual and Jacob...
double side_length(unsigned int side) const
void add_contribution(const T *local, T **y_global) const
Add the values of element-local contributions y to the global vector y_global.
void nodal_values(const array::Scalar &x_global, int *result) const
Get nodal values of an integer mask.
unsigned int n_sides() const
Vector2d normal(unsigned int side) const
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.
double weight(unsigned int q) const
Weight of the quadrature point q
const Germ & chi(unsigned int q, unsigned int k) const
unsigned int n_pts() const
Number of quadrature points.
The mapping from global to local degrees of freedom.
P1 element embedded in Q1Element2.
double weight(int k) const
QuadPoint point(int k) const
std::shared_ptr< TerminationReason > solve_with_reason(const Inputs &inputs)
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 ...
void driving_stress(const fem::Element &E, const Coefficients *x, Vector2d *result) const
const array::Scalar * m_driving_stress_x
void cache_inputs(const Inputs &inputs)
Initialize stored data from the coefficients in the SSA. Called by SSAFEM::solve.
static PetscErrorCode function_callback(DMDALocalInfo *info, Vector2d const *const *velocity, Vector2d **residual, CallbackData *fe)
SNES callbacks.
double m_beta_ice_free_bedrock
SSAFEM(std::shared_ptr< const Grid > g)
static void explicit_driving_stress(const fem::Element &E, const Coefficients *x, Vector2d *result)
CallbackData m_callback_data
void compute_local_jacobian(Vector2d const *const *velocity, Mat J)
Implements the callback for computing the Jacobian.
array::Vector1 m_boundary_integral
Boundary integral (CFBC contribution to the residual).
static PetscErrorCode jacobian_callback(DMDALocalInfo *info, Vector2d const *const *velocity, Mat A, Mat J, CallbackData *fe)
void monitor_jacobian(Mat Jac)
void monitor_function(Vector2d const *const *velocity_global, Vector2d const *const *residual_global)
virtual void solve(const Inputs &inputs)
fem::Q1Element2 m_q1_element
virtual void init_impl()
Initialize a generic regular-grid SSA solver.
std::shared_ptr< TerminationReason > solve_nocache()
void cache_residual_cfbc(const Inputs &inputs)
Compute and cache residual contributions from the integral over the lateral boundary.
void compute_local_function(Vector2d const *const *velocity, Vector2d **residual)
Implements the callback for computing the residual.
fem::ElementIterator m_element_index
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.
const array::Scalar * m_driving_stress_y
array::Scalar1 m_node_type
Storage for node types (interior, boundary, exterior).
array::Vector1 m_bc_values
array::Array2D< Coefficients > m_coefficients
double get_notional_strength() const
Returns strength = (viscosity times thickness).
array::Vector m_velocity_global
SSAStrengthExtension * strength_extension
virtual void init_impl()
Initialize a generic regular-grid SSA solver.
const array::Vector1 & velocity() const
Get the thickness-advective 2D velocity.
std::shared_ptr< rheology::FlowLaw > m_flow_law
IceBasalResistancePlasticLaw * m_basal_sliding_law
array::Vector2 m_velocity
#define PISM_CHK(errcode, name)
#define PISM_ERROR_LOCATION
Germ chi(unsigned int k, const QuadPoint &pt)
Linear basis functions on the interval [-1, -1].
const unsigned int MAX_QUADRATURE_SIZE
ElementType element_type(const int node_type[q1::n_chi])
bool ice_free_land(int M)
bool grounded(int M)
Grounded cell (grounded ice or ice-free).
bool ocean(int M)
An ocean cell (floating ice or ice-free).
bool Bool(const std::string &option, const std::string &description)
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 .
static double secondInvariant_2D(const Vector2d &U_x, const Vector2d &U_y)
void compute_node_types(const array::Scalar1 &ice_thickness, double thickness_threshold, array::Scalar &result)
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,...)
void handle_fatal_errors(MPI_Comm com)
double dy
Function derivative with respect to y.
double val
Function value.
double dx
Function deriviative with respect to x.
Struct for gathering the value and derivative of a function at a point.
Adaptor for gluing SNESDAFormFunction callbacks to an SSAFEM.
double thickness
ice thickness
double tauc
basal yield stress
double sea_level
sea level
Vector2d driving_stress
prescribed gravitational driving stress
double hardness
ice hardness