19#include "pism/inverse/IP_SSAHardavForwardProblem.hh"
20#include "pism/util/Grid.hh"
21#include "pism/util/Config.hh"
22#include "pism/util/Vars.hh"
23#include "pism/util/error_handling.hh"
24#include "pism/rheology/FlowLaw.hh"
25#include "pism/geometry/Geometry.hh"
26#include "pism/stressbalance/StressBalance.hh"
27#include "pism/util/petscwrappers/DM.hh"
28#include "pism/util/petscwrappers/Vec.hh"
29#include "pism/util/fem/Quadrature.hh"
30#include "pism/util/fem/DirichletData.hh"
31#include "pism/util/Logger.hh"
41 m_dzeta_local(m_grid,
"d_zeta_local"),
42 m_fixed_design_locations(NULL),
44 m_du_global(m_grid,
"linearization work vector (sans ghosts)"),
45 m_du_local(m_grid,
"linearization work vector (with ghosts)"),
46 m_hardav(m_grid,
"hardav"),
47 m_element_index(*m_grid),
48 m_element(*m_grid, fem::Q1Quadrature4()),
49 m_rebuild_J_state(true)
60 ierr = DMSetMatType(*dm, MATBAIJ);
69 double ksp_rtol = 1e-12;
70 ierr = KSPSetTolerances(
m_ksp, ksp_rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
74 ierr = KSPGetPC(
m_ksp, &pc);
77 ierr = PCSetType(pc, PCBJACOBI);
80 ierr = KSPSetFromOptions(
m_ksp);
97 if (
m_config->get_flag(
"geometry.part_grid.enabled")) {
99 *
m_grid->variables().get_2d_scalar(
"ice_area_specific_volume"));
108 const auto &variables =
m_grid->variables();
111 if (variables.is_available(
"vel_bc_mask")) {
112 vel_bc_mask = variables.get_2d_scalar(
"vel_bc_mask");
116 if (variables.is_available(
"vel_bc")) {
117 vel_bc = variables.get_2d_vector(
"vel_bc");
123 inputs.
enthalpy = variables.get_3d_scalar(
"enthalpy");
153 for (
auto p :
m_grid->points_with_ghosts(1)) {
154 const int i = p.i(), j = p.j();
262 dzeta_local = &dzeta;
267 list.
add(*dzeta_local);
270 for (
auto p :
m_grid->points()) {
271 const int i = p.i(), j = p.j();
283 Vector2d U[Nq_max], U_x[Nq_max], U_y[Nq_max];
307 for (
int j =ys; j<ys+ym; j++) {
308 for (
int i =xs; i<xs+xm; i++) {
311 for (
unsigned int k=0;
k<Nk;
k++) {
336 for (
unsigned int k=0;
k<Nk;
k++) {
342 double thickness[Nq_max];
347 double hardness[Nq_max];
352 mask, thickness, tauc, hardness);
355 for (
unsigned int q = 0; q < Nq; q++) {
357 double Duqq[3] = {U_x[q].
u, U_y[q].
v, 0.5 * (U_y[q].
u + U_x[q].
v)};
364 d_nuH *= (2.0 * thickness[q]);
369 for (
unsigned int k = 0;
k < Nk;
k++) {
371 du_e[
k].
u += W*d_nuH*(testqk.
dx*(2*Duqq[0] + Duqq[1]) + testqk.
dy*Duqq[2]);
372 du_e[
k].
v += W*d_nuH*(testqk.
dy*(2*Duqq[1] + Duqq[0]) + testqk.
dx*Duqq[2]);
406 auto da2 =
m_grid->get_dm(1,
m_config->get_number(
"grid.max_stencil_width"));
451 Vector2d U[Nq_max], U_x[Nq_max], U_y[Nq_max];
469 for (
auto p :
m_grid->points()) {
470 const int i = p.i(), j = p.j();
483 for (
int j = ys; j < ys + ym; j++) {
484 for (
int i = xs; i < xs + xm; i++) {
503 for (
unsigned int k = 0;
k < Nk;
k++) {
507 double thickness[Nq_max];
512 double hardness[Nq_max];
517 mask, thickness, tauc, hardness);
520 for (
unsigned int q = 0; q < Nq; q++) {
522 double Duqq[3] = {U_x[q].
u, U_y[q].
v, 0.5 * (U_y[q].
u + U_x[q].
v)};
530 d_nuH_dB *= (2.0 * thickness[q]);
535 for (
unsigned int k = 0;
k < Nk;
k++) {
536 dzeta_e[
k] += W*d_nuH_dB*
m_element.
chi(q,
k).
val*((du_dx_q[q].
u*(2*Duqq[0] + Duqq[1]) +
537 du_dy_q[q].u*Duqq[2]) +
538 (du_dy_q[q].
v*(2*Duqq[1] + Duqq[0]) +
539 du_dx_q[q].v*Duqq[2]));
551 for (
auto p :
m_grid->points()) {
552 const int i = p.i(), j = p.j();
556 dzeta_a[j][i] *= dB_dzeta;
597 KSPConvergedReason reason;
598 ierr = KSPGetConvergedReason(
m_ksp, &reason);
599 PISM_CHK(ierr,
"KSPGetConvergedReason");
603 " failed to converge (KSP reason %s)",
604 KSPConvergedReasons[reason]);
608 "IP_SSAHardavForwardProblem::apply_linearization converged"
609 " (KSP reason %s)\n",
610 KSPConvergedReasons[reason]);
663 KSPConvergedReason reason;
664 ierr = KSPGetConvergedReason(
m_ksp, &reason);
665 PISM_CHK(ierr,
"KSPGetConvergedReason");
669 KSPConvergedReasons[reason]);
672 m_log->message(4,
"IP_SSAHardavForwardProblem::apply_linearization converged (KSP reason %s)\n",
673 KSPConvergedReasons[reason]);
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)
array::Scalar1 sea_level_elevation
void ensure_consistency(double ice_free_thickness_threshold)
array::Scalar1 ice_area_specific_volume
array::Scalar2 ice_thickness
array::Scalar2 bed_elevation
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)
void scale(double alpha)
Result: v <- v * alpha. Calls VecScale.
void set(double c)
Result: v[j] <- c for all j.
std::shared_ptr< petsc::DM > dm() const
void update_ghosts()
Updates ghost points.
VariableMetadata & metadata(unsigned int N=0)
Returns a reference to the VariableMetadata object containing metadata for the compoment N.
unsigned int stencil_width() const
Get the stencil width of the current Array. Returns 0 if ghosts are not available.
void enforce_homogeneous(const Element2 &element, double *x_e)
void fix_residual_homogeneous(double **r_global)
void enforce(const Element2 &element, Vector2d *x_e)
void fix_residual_homogeneous(Vector2d **r)
void enforce_homogeneous(const Element2 &element, Vector2d *x_e)
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...
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.
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.
virtual void convertToDesignVariable(array::Scalar &zeta, array::Scalar &d, bool communicate=true)
Transforms a vector of values to a vector of values.
virtual void toDesignVariable(double zeta, double *value, double *derivative)=0
Converts from parameterization value to .
petsc::Mat m_J_state
Mat used in apply_linearization and apply_linearization_transpose.
virtual void assemble_residual(array::Vector &u, array::Vector &R)
Computes the residual function as defined in the class-level documentation.
virtual void apply_jacobian_design_transpose(array::Vector &u, array::Vector &du, array::Scalar &dzeta)
Applies the transpose of the design Jacobian matrix to a perturbation of the state variable.
array::Vector m_du_global
Temporary storage when state vectors need to be used without ghosts.
virtual void apply_linearization(array::Scalar &dzeta, array::Vector &du)
Applies the linearization of the forward map (i.e. the reduced gradient described in the class-level...
virtual void assemble_jacobian_state(array::Vector &u, Mat J)
Assembles the state Jacobian matrix.
IPDesignVariableParameterization & m_design_param
The function taking to .
array::Scalar1 m_hardav
Vertically-averaged ice hardness.
array::Scalar * m_fixed_design_locations
Locations where should not be adjusted.
array::Scalar1 m_dzeta_local
Storage for d_zeta with ghosts, if needed when an argument d_zeta is ghost-less.
IP_SSAHardavForwardProblem(std::shared_ptr< const Grid > g, IPDesignVariableParameterization &tp)
Constructs from the same objects as SSAFEM, plus a specification of how is parameterized.
fem::ElementIterator m_element_index
std::shared_ptr< array::Vector > m_velocity_shared
virtual void apply_jacobian_design(array::Vector &u, array::Scalar &dzeta, array::Vector &du)
Applies the design Jacobian matrix to a perturbation of the design variable.
virtual std::shared_ptr< TerminationReason > linearize_at(array::Scalar &zeta)
Sets the current value of the design variable and solves the SSA to find the associated .
array::Scalar * m_zeta
Current value of zeta, provided from caller.
bool m_rebuild_J_state
Flag indicating that the state jacobian matrix needs rebuilding.
virtual void apply_linearization_transpose(array::Vector &du, array::Scalar &dzeta)
Applies the transpose of the linearization of the forward map (i.e. the transpose of the reduced grad...
virtual void set_design(array::Scalar &zeta)
Sets the current value of of the design paramter .
fem::Q1Element2 m_element
array::Vector1 m_du_local
Temporary storage when state vectors need to be used with ghosts.
petsc::KSP m_ksp
KSP used in apply_linearization and apply_linearization_transpose.
void cache_inputs(const Inputs &inputs)
Initialize stored data from the coefficients in the SSA. Called by SSAFEM::solve.
void compute_local_jacobian(Vector2d const *const *velocity, Mat J)
Implements the callback for computing the Jacobian.
std::shared_ptr< TerminationReason > solve_nocache()
void compute_local_function(Vector2d const *const *velocity, Vector2d **residual)
Implements the callback for computing the residual.
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.
array::Vector1 m_bc_values
array::Array2D< Coefficients > m_coefficients
double get_min_thickness() const
Returns minimum thickness to trigger use of extension.
array::Vector m_velocity_global
SSAStrengthExtension * strength_extension
std::shared_ptr< rheology::FlowLaw > m_flow_law
array::Vector2 m_velocity
#define PISM_CHK(errcode, name)
#define PISM_ERROR_LOCATION
const unsigned int MAX_QUADRATURE_SIZE
static double secondInvariant_2D(const Vector2d &U_x, const Vector2d &U_y)
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.