PISM, A Parallel Ice Sheet Model 2.3.0-79cae578d committed by Constantine Khrulev on 2026-03-22
Loading...
Searching...
No Matches
IP_SSAHardavForwardProblem.cc
Go to the documentation of this file.
1// Copyright (C) 2013, 2014, 2015, 2016, 2017, 2018, 2020, 2021, 2022, 2023, 2024, 2025 David Maxwell and Constantine Khroulev
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#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"
32
33namespace pism {
34namespace inverse {
35
38 : SSAFEM(g),
39 m_stencil_width(1),
40 m_zeta(NULL),
41 m_dzeta_local(m_grid, "d_zeta_local"),
42 m_fixed_design_locations(NULL),
43 m_design_param(tp),
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)
50{
51
52 PetscErrorCode ierr;
53
54 m_velocity_shared.reset(new array::Vector(m_grid, "dummy"));
55 m_velocity_shared->metadata(0) = m_velocity.metadata(0);
56 m_velocity_shared->metadata(1) = m_velocity.metadata(1);
57
58 auto dm = m_velocity_global.dm();
59
60 ierr = DMSetMatType(*dm, MATBAIJ);
61 PISM_CHK(ierr, "DMSetMatType");
62
63 ierr = DMCreateMatrix(*dm, m_J_state.rawptr());
64 PISM_CHK(ierr, "DMCreateMatrix");
65
66 ierr = KSPCreate(m_grid->com, m_ksp.rawptr());
67 PISM_CHK(ierr, "KSPCreate");
68
69 double ksp_rtol = 1e-12;
70 ierr = KSPSetTolerances(m_ksp, ksp_rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
71 PISM_CHK(ierr, "KSPSetTolerances");
72
73 PC pc;
74 ierr = KSPGetPC(m_ksp, &pc);
75 PISM_CHK(ierr, "KSPGetPC");
76
77 ierr = PCSetType(pc, PCBJACOBI);
78 PISM_CHK(ierr, "PCSetType");
79
80 ierr = KSPSetFromOptions(m_ksp);
81 PISM_CHK(ierr, "KSPSetFromOptions");
82}
83
85
86 SSAFEM::init();
87
88 // Get most of the inputs from Grid::variables() and fake the rest.
89 //
90 // I will need to fix this at some point.
91 {
92 Geometry geometry(m_grid);
93 geometry.ice_thickness.copy_from(*m_grid->variables().get_2d_scalar("land_ice_thickness"));
94 geometry.bed_elevation.copy_from(*m_grid->variables().get_2d_scalar("bedrock_altitude"));
95 geometry.sea_level_elevation.set(0.0); // FIXME: this should be an input
96
97 if (m_config->get_flag("geometry.part_grid.enabled")) {
99 *m_grid->variables().get_2d_scalar("ice_area_specific_volume"));
100 } else {
101 geometry.ice_area_specific_volume.set(0.0);
102 }
103
104 geometry.ensure_consistency(m_config->get_number("stress_balance.ice_free_thickness_standard"));
105
107
108 const auto &variables = m_grid->variables();
109
110 const array::Scalar *vel_bc_mask = nullptr;
111 if (variables.is_available("vel_bc_mask")) {
112 vel_bc_mask = variables.get_2d_scalar("vel_bc_mask");
113 }
114
115 const array::Vector *vel_bc = nullptr;
116 if (variables.is_available("vel_bc")) {
117 vel_bc = variables.get_2d_vector("vel_bc");
118 }
119
120 inputs.geometry = &geometry;
121 inputs.basal_melt_rate = NULL;
122 inputs.basal_yield_stress = variables.get_2d_scalar("tauc");
123 inputs.enthalpy = variables.get_3d_scalar("enthalpy");
124 inputs.age = NULL;
125 inputs.bc_mask = vel_bc_mask;
126 inputs.bc_values = vel_bc;
127
128 inputs.water_column_pressure = NULL;
129
130 cache_inputs(inputs);
131 }
132}
133
134//! Sets the current value of of the design paramter \f$\zeta\f$.
135/*! This method sets \f$\zeta\f$ but does not solve the %SSA.
136It it intended for inverse methods that simultaneously compute
137the pair \f$u\f$ and \f$\zeta\f$ without ever solving the %SSA
138directly. Use this method in conjuction with
139\ref assemble_jacobian_state and \ref apply_jacobian_design and their friends.
140The vector \f$\zeta\f$ is not copied; a reference to the array::Array is
141kept.
142*/
144
145 m_zeta = &new_zeta;
146
147 // Convert zeta to hardav.
149
150 // Cache hardav at the quadrature points.
152
153 for (auto p : m_grid->points_with_ghosts(1)) {
154 const int i = p.i(), j = p.j();
155 m_coefficients(i, j).hardness = m_hardav(i, j);
156 }
157
158 // Flag the state jacobian as needing rebuilding.
159 m_rebuild_J_state = true;
160}
161
162//! Sets the current value of the design variable \f$\zeta\f$ and solves the %SSA to find the associated \f$u_{\rm SSA}\f$.
163/* Use this method for inverse methods employing the reduced gradient. Use this method
164in conjuction with apply_linearization and apply_linearization_transpose.*/
165std::shared_ptr<TerminationReason> IP_SSAHardavForwardProblem::linearize_at(array::Scalar &zeta) {
166 this->set_design(zeta);
167 return this->solve_nocache();
168}
169
170//! Computes the residual function \f$\mathcal{R}(u, \zeta)\f$ as defined in the class-level documentation.
171/* The value of \f$\zeta\f$ is set prior to this call via set_design or linearize_at. The value
172of the residual is returned in \a RHS.*/
179
180//! Computes the residual function \f$\mathcal{R}(u, \zeta)\f$ defined in the class-level documentation.
181/* The return value is specified via a Vec for the benefit of certain TAO routines. Otherwise,
182the method is identical to the assemble_residual returning values as a StateVec (an array::Vector).*/
184
185 array::AccessScope l{&u};
186
188
189 this->compute_local_function(u.array(), (Vector2d**)rhs_a.get());
190}
191
192//! Assembles the state Jacobian matrix.
193/* The matrix depends on the current value of the design variable \f$\zeta\f$ and the current
194value of the state variable \f$u\f$. The specification of \f$\zeta\f$ is done earlier
195with set_design or linearize_at. The value of \f$u\f$ is specified explicitly as an argument
196to this method.
197 @param[in] u Current state variable value.
198 @param[out] J computed state Jacobian.
199*/
205
206//! Applies the design Jacobian matrix to a perturbation of the design variable.
207/*! The return value uses a DesignVector (array::Vector), which can be ghostless. Ghosts (if present) are updated.
208\overload
209*/
217
218//! Applies the design Jacobian matrix to a perturbation of the design variable.
219/*! The return value is a Vec for the benefit of TAO. It is assumed to be ghostless; no communication is done.
220\overload
221*/
228
229//! @brief Applies the design Jacobian matrix to a perturbation of the
230//! design variable.
231
232/*! The matrix depends on the current value of the design variable
233 \f$\zeta\f$ and the current value of the state variable \f$u\f$.
234 The specification of \f$\zeta\f$ is done earlier with set_design
235 or linearize_at. The value of \f$u\f$ is specified explicitly as
236 an argument to this method.
237
238 @param[in] u Current state variable value.
239
240 @param[in] dzeta Perturbation of the design variable. Prefers
241 vectors with ghosts; will copy to a ghosted vector
242 if needed.
243
244 @param[out] du_a Computed corresponding perturbation of the state
245 variable. The array \a du_a should be extracted
246 first from a Vec or an array::Array.
247
248 Typically this method is called via one of its overloads.
249*/
251 array::Scalar &dzeta,
252 Vector2d **du_a) {
253
254 const unsigned int Nk = fem::q1::n_chi;
255 const unsigned int Nq = m_element.n_pts();
256 const unsigned int Nq_max = fem::MAX_QUADRATURE_SIZE;
257
259
260 array::Scalar *dzeta_local;
261 if (dzeta.stencil_width() > 0) {
262 dzeta_local = &dzeta;
263 } else {
265 dzeta_local = &m_dzeta_local;
266 }
267 list.add(*dzeta_local);
268
269 // Zero out the portion of the function we are responsible for computing.
270 for (auto p : m_grid->points()) {
271 const int i = p.i(), j = p.j();
272
273 du_a[j][i].u = 0.0;
274 du_a[j][i].v = 0.0;
275 }
276
277 // Aliases to help with notation consistency below.
278 const array::Scalar *dirichletLocations = &m_bc_mask;
279 const array::Vector *dirichletValues = &m_bc_values;
280 double dirichletWeight = m_dirichletScale;
281
282 Vector2d u_e[Nk];
283 Vector2d U[Nq_max], U_x[Nq_max], U_y[Nq_max];
284
285 Vector2d du_e[Nk];
286
287 double dzeta_e[Nk];
288
289 double zeta_e[Nk];
290
291 double dB_e[Nk];
292 double dB_q[Nq_max];
293
294 fem::DirichletData_Vector dirichletBC(dirichletLocations, dirichletValues,
295 dirichletWeight);
297
298 // Loop through all elements.
299 const int
300 xs = m_element_index.xs,
301 xm = m_element_index.xm,
302 ys = m_element_index.ys,
303 ym = m_element_index.ym;
304
305 ParallelSection loop(m_grid->com);
306 try {
307 for (int j =ys; j<ys+ym; j++) {
308 for (int i =xs; i<xs+xm; i++) {
309
310 // Zero out the element-local residual in prep for updating it.
311 for (unsigned int k=0; k<Nk; k++) {
312 du_e[k].u = 0;
313 du_e[k].v = 0;
314 }
315
316 // Initialize the map from global to local degrees of freedom for this element.
317 m_element.reset(i, j);
318
319 // Obtain the value of the solution at the nodes adjacent to the element,
320 // fix dirichlet values, and compute values at quad pts.
321 m_element.nodal_values(u.array(), u_e);
322 if (dirichletBC) {
323 dirichletBC.constrain(m_element);
324 dirichletBC.enforce(m_element, u_e);
325 }
326 m_element.evaluate(u_e, U, U_x, U_y);
327
328 // Compute dzeta at the nodes
329 m_element.nodal_values(dzeta_local->array(), dzeta_e);
330 if (fixedZeta) {
331 fixedZeta.enforce_homogeneous(m_element, dzeta_e);
332 }
333
334 // Compute the change in hardav with respect to zeta at the quad points.
336 for (unsigned int k=0; k<Nk; k++) {
337 m_design_param.toDesignVariable(zeta_e[k], NULL, dB_e + k);
338 dB_e[k]*=dzeta_e[k];
339 }
340 m_element.evaluate(dB_e, dB_q);
341
342 double thickness[Nq_max];
343 {
344 Coefficients coeffs[Nk];
345 int mask[Nq_max];
346 double tauc[Nq_max];
347 double hardness[Nq_max];
348
349 m_element.nodal_values(m_coefficients.array(), coeffs);
350
352 mask, thickness, tauc, hardness);
353 }
354
355 for (unsigned int q = 0; q < Nq; q++) {
356 // Symmetric gradient at the quadrature point.
357 double Duqq[3] = {U_x[q].u, U_y[q].v, 0.5 * (U_y[q].u + U_x[q].v)};
358
359 double d_nuH = 0;
360 if (thickness[q] >= strength_extension->get_min_thickness()) {
361 m_flow_law->effective_viscosity(dB_q[q],
362 secondInvariant_2D(U_x[q], U_y[q]),
363 &d_nuH, NULL);
364 d_nuH *= (2.0 * thickness[q]);
365 }
366
367 auto W = m_element.weight(q);
368
369 for (unsigned int k = 0; k < Nk; k++) {
370 const fem::Germ &testqk = m_element.chi(q, 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]);
373 }
374 } // q
375 m_element.add_contribution(du_e, du_a);
376 } // j
377 } // i
378 } catch (...) {
379 loop.failed();
380 }
381 loop.check();
382
383 if (dirichletBC) {
384 dirichletBC.fix_residual_homogeneous(du_a);
385 }
386}
387
388//! Applies the transpose of the design Jacobian matrix to a perturbation of the state variable.
389/*! The return value uses a StateVector (array::Scalar) which can be ghostless; ghosts (if present) are updated.
390\overload
391*/
398
399//! Applies the transpose of the design Jacobian matrix to a perturbation of the state variable.
400/*! The return value uses a Vec for the benefit of TAO. It is assumed to be ghostless; no communication is done.
401\overload */
403 array::Vector &du,
404 Vec dzeta) {
405
406 auto da2 = m_grid->get_dm(1, m_config->get_number("grid.max_stencil_width"));
407 petsc::DMDAVecArray dzeta_a(da2, dzeta);
408 this->apply_jacobian_design_transpose(u, du, (double**)dzeta_a.get());
409}
410
411//! @brief Applies the transpose of the design Jacobian matrix to a
412//! perturbation of the state variable.
413
414/*! The matrix depends on the current value of the design variable
415 \f$\zeta\f$ and the current value of the state variable \f$u\f$.
416 The specification of \f$\zeta\f$ is done earlier with set_design
417 or linearize_at. The value of \f$u\f$ is specified explicitly as
418 an argument to this method.
419
420 @param[in] u Current state variable value.
421
422 @param[in] du Perturbation of the state variable. Prefers vectors
423 with ghosts; will copy to a ghosted vector if need be.
424
425 @param[out] dzeta_a Computed corresponding perturbation of the
426 design variable. The array \a dzeta_a should be
427 extracted first from a Vec or an array::Array.
428
429 Typically this method is called via one of its overloads.
430*/
432 array::Vector &du,
433 double **dzeta_a) {
434
435 const unsigned int Nk = fem::q1::n_chi;
436 const unsigned int Nq = m_element.n_pts();
437 const unsigned int Nq_max = fem::MAX_QUADRATURE_SIZE;
438
440
441 array::Vector *du_local;
442 if (du.stencil_width() > 0) {
443 du_local = &du;
444 } else {
446 du_local = &m_du_local;
447 }
448 list.add(*du_local);
449
450 Vector2d u_e[Nk];
451 Vector2d U[Nq_max], U_x[Nq_max], U_y[Nq_max];
452
453 Vector2d du_e[Nk];
454 Vector2d du_q[Nq_max];
455 Vector2d du_dx_q[Nq_max];
456 Vector2d du_dy_q[Nq_max];
457
458 double dzeta_e[Nk];
459
460 // Aliases to help with notation consistency.
461 const array::Scalar *dirichletLocations = &m_bc_mask;
462 const array::Vector *dirichletValues = &m_bc_values;
463 double dirichletWeight = m_dirichletScale;
464
465 fem::DirichletData_Vector dirichletBC(dirichletLocations, dirichletValues,
466 dirichletWeight);
467
468 // Zero out the portion of the function we are responsible for computing.
469 for (auto p : m_grid->points()) {
470 const int i = p.i(), j = p.j();
471
472 dzeta_a[j][i] = 0;
473 }
474
475 const int
476 xs = m_element_index.xs,
477 xm = m_element_index.xm,
478 ys = m_element_index.ys,
479 ym = m_element_index.ym;
480
481 ParallelSection loop(m_grid->com);
482 try {
483 for (int j = ys; j < ys + ym; j++) {
484 for (int i = xs; i < xs + xm; i++) {
485 // Initialize the map from global to local degrees of freedom for this element.
486 m_element.reset(i, j);
487
488 // Obtain the value of the solution at the nodes adjacent to the element.
489 // Compute the solution values and symmetric gradient at the quadrature points.
490 m_element.nodal_values(du.array(), du_e);
491 if (dirichletBC) {
492 dirichletBC.enforce_homogeneous(m_element, du_e);
493 }
494 m_element.evaluate(du_e, du_q, du_dx_q, du_dy_q);
495
496 m_element.nodal_values(u.array(), u_e);
497 if (dirichletBC) {
498 dirichletBC.enforce(m_element, u_e);
499 }
500 m_element.evaluate(u_e, U, U_x, U_y);
501
502 // Zero out the element-local residual in prep for updating it.
503 for (unsigned int k = 0; k < Nk; k++) {
504 dzeta_e[k] = 0;
505 }
506
507 double thickness[Nq_max];
508 {
509 Coefficients coeffs[Nk];
510 int mask[Nq_max];
511 double tauc[Nq_max];
512 double hardness[Nq_max];
513
514 m_element.nodal_values(m_coefficients.array(), coeffs);
515
517 mask, thickness, tauc, hardness);
518 }
519
520 for (unsigned int q = 0; q < Nq; q++) {
521 // Symmetric gradient at the quadrature point.
522 double Duqq[3] = {U_x[q].u, U_y[q].v, 0.5 * (U_y[q].u + U_x[q].v)};
523
524 // Determine "d_nuH / dB" at the quadrature point
525 double d_nuH_dB = 0;
526 if (thickness[q] >= strength_extension->get_min_thickness()) {
527 m_flow_law->effective_viscosity(1.0,
528 secondInvariant_2D(U_x[q], U_y[q]),
529 &d_nuH_dB, NULL);
530 d_nuH_dB *= (2.0 * thickness[q]);
531 }
532
533 auto W = m_element.weight(q);
534
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]));
540 }
541 } // q
542
543 m_element.add_contribution(dzeta_e, dzeta_a);
544 } // j
545 } // i
546 } catch (...) {
547 loop.failed();
548 }
549 loop.check();
550
551 for (auto p : m_grid->points()) {
552 const int i = p.i(), j = p.j();
553
554 double dB_dzeta;
555 m_design_param.toDesignVariable((*m_zeta)(i, j), NULL, &dB_dzeta);
556 dzeta_a[j][i] *= dB_dzeta;
557 }
558
559 if (m_fixed_design_locations != nullptr) {
561 fixedZeta.fix_residual_homogeneous(dzeta_a);
562 }
563}
564
565/*!\brief Applies the linearization of the forward map (i.e. the reduced gradient \f$DF\f$ described in
566the class-level documentation.) */
567/*! As described previously,
568\f[
569Df = J_{\rm State}^{-1} J_{\rm Design}.
570\f]
571Applying the linearization then involves the solution of a linear equation.
572The matrices \f$J_{\rm State}\f$ and \f$J_{\rm Design}\f$ both depend on the value of the
573design variable \f$\zeta\f$ and the value of the corresponding state variable \f$u=F(\zeta)\f$.
574These are established by first calling linearize_at.
575 @param[in] dzeta Perturbation of the design variable
576 @param[out] du Computed corresponding perturbation of the state variable; ghosts (if present) are updated.
577*/
579
580 PetscErrorCode ierr;
581
582 if (m_rebuild_J_state) {
584 m_rebuild_J_state = false;
585 }
586
588 m_du_global.scale(-1);
589
590 // call PETSc to solve linear system by iterative method.
591 ierr = KSPSetOperators(m_ksp, m_J_state, m_J_state);
592 PISM_CHK(ierr, "KSPSetOperators");
593
594 ierr = KSPSolve(m_ksp, m_du_global.vec(), m_du_global.vec());
595 PISM_CHK(ierr, "KSPSolve"); // SOLVE
596
597 KSPConvergedReason reason;
598 ierr = KSPGetConvergedReason(m_ksp, &reason);
599 PISM_CHK(ierr, "KSPGetConvergedReason");
600
601 if (reason < 0) {
602 throw RuntimeError::formatted(PISM_ERROR_LOCATION, "IP_SSAHardavForwardProblem::apply_linearization solve"
603 " failed to converge (KSP reason %s)",
604 KSPConvergedReasons[reason]);
605 }
606
607 m_log->message(4,
608 "IP_SSAHardavForwardProblem::apply_linearization converged"
609 " (KSP reason %s)\n",
610 KSPConvergedReasons[reason]);
611
613}
614
615/*! \brief Applies the transpose of the linearization of the forward map
616 (i.e. the transpose of the reduced gradient \f$DF\f$ described in the class-level documentation.) */
617/*! As described previously,
618\f[
619Df = J_{\rm State}^{-1} J_{\rm Design}.
620\f]
621so
622\f[
623Df^t = J_{\rm Design}^t \; (J_{\rm State}^t)^{-1} .
624\f]
625Applying the transpose of the linearization then involves the solution of a linear equation.
626The matrices \f$J_{\rm State}\f$ and \f$J_{\rm Design}\f$ both depend on the value of the
627design variable \f$\zeta\f$ and the value of the corresponding state variable \f$u=F(\zeta)\f$.
628These are established by first calling linearize_at.
629 @param[in] du Perturbation of the state variable
630 @param[out] dzeta Computed corresponding perturbation of the design variable; ghosts (if present) are updated.
631*/
633 array::Scalar &dzeta) {
634
635 PetscErrorCode ierr;
636
637 if (m_rebuild_J_state) {
639 m_rebuild_J_state = false;
640 }
641
642 // Aliases to help with notation consistency below.
643 const array::Scalar *dirichletLocations = &m_bc_mask;
644 const array::Vector *dirichletValues = &m_bc_values;
645 double dirichletWeight = m_dirichletScale;
646
648
649 fem::DirichletData_Vector dirichletBC(dirichletLocations, dirichletValues,
650 dirichletWeight);
651 if (dirichletBC) {
654 }
655
656 // call PETSc to solve linear system by iterative method.
657 ierr = KSPSetOperators(m_ksp, m_J_state, m_J_state);
658 PISM_CHK(ierr, "KSPSetOperators");
659
660 ierr = KSPSolve(m_ksp, m_du_global.vec(), m_du_global.vec());
661 PISM_CHK(ierr, "KSPSolve"); // SOLVE
662
663 KSPConvergedReason reason;
664 ierr = KSPGetConvergedReason(m_ksp, &reason);
665 PISM_CHK(ierr, "KSPGetConvergedReason");
666
667 if (reason < 0) {
668 throw RuntimeError::formatted(PISM_ERROR_LOCATION, "IP_SSAHardavForwardProblem::apply_linearization solve failed to converge (KSP reason %s)",
669 KSPConvergedReasons[reason]);
670 }
671
672 m_log->message(4, "IP_SSAHardavForwardProblem::apply_linearization converged (KSP reason %s)\n",
673 KSPConvergedReasons[reason]);
674
676 dzeta.scale(-1);
677
678 if (dzeta.stencil_width() > 0) {
679 dzeta.update_ghosts();
680 }
681}
682
683} // end of namespace inverse
684} // end of namespace pism
std::shared_ptr< const Config > m_config
configuration database used by this component
Definition Component.hh:160
const std::shared_ptr< const Grid > m_grid
grid used by this component
Definition Component.hh:158
std::shared_ptr< const Logger > m_log
logger (for easy access)
Definition Component.hh:164
array::Scalar1 sea_level_elevation
Definition Geometry.hh:48
void ensure_consistency(double ice_free_thickness_threshold)
Definition Geometry.cc:116
array::Scalar1 ice_area_specific_volume
Definition Geometry.hh:52
array::Scalar2 ice_thickness
Definition Geometry.hh:51
array::Scalar2 bed_elevation
Definition Geometry.hh:47
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.
Definition Vector2d.hh:29
T * rawptr()
Definition Wrapper.hh:34
Makes sure that we call begin_access() and end_access() for all accessed array::Arrays.
Definition Array.hh:66
void copy_from(const Array2D< T > &source)
Definition Array2D.hh:101
void add(double alpha, const Array2D< T > &x)
Definition Array2D.hh:93
petsc::Vec & vec() const
Definition Array.cc:313
void scale(double alpha)
Result: v <- v * alpha. Calls VecScale.
Definition Array.cc:227
void set(double c)
Result: v[j] <- c for all j.
Definition Array.cc:659
std::shared_ptr< petsc::DM > dm() const
Definition Array.cc:327
void update_ghosts()
Updates ghost points.
Definition Array.cc:645
VariableMetadata & metadata(unsigned int N=0)
Returns a reference to the VariableMetadata object containing metadata for the compoment N.
Definition Array.cc:467
unsigned int stencil_width() const
Get the stencil width of the current Array. Returns 0 if ghosts are not available.
Definition Array.cc:305
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.
Definition Element.hh:195
void reset(int i, int j)
Initialize the Element to element (i, j) for the purposes of inserting into global residual and Jacob...
Definition Element.cc:196
void add_contribution(const T *local, T **y_global) const
Add the values of element-local contributions y to the global vector y_global.
Definition Element.hh:238
void nodal_values(const array::Scalar &x_global, int *result) const
Get nodal values of an integer mask.
Definition Element.cc:185
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
Definition Element.hh:85
const Germ & chi(unsigned int q, unsigned int k) const
Definition Element.hh:73
unsigned int n_pts() const
Number of quadrature points.
Definition Element.hh:80
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.
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 .
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.
const array::Scalar * basal_yield_stress
const array::Scalar * bc_mask
const array::Array3D * age
const array::Scalar * water_column_pressure
const array::Array3D * enthalpy
const array::Scalar * basal_melt_rate
const array::Vector * bc_values
void cache_inputs(const Inputs &inputs)
Initialize stored data from the coefficients in the SSA. Called by SSAFEM::solve.
Definition SSAFEM.cc:268
void compute_local_jacobian(Vector2d const *const *velocity, Mat J)
Implements the callback for computing the Jacobian.
Definition SSAFEM.cc:958
array::Scalar1 m_bc_mask
Definition SSAFEM.hh:67
std::shared_ptr< TerminationReason > solve_nocache()
Definition SSAFEM.cc:195
void compute_local_function(Vector2d const *const *velocity, Vector2d **residual)
Implements the callback for computing the residual.
Definition SSAFEM.cc:736
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:346
array::Vector1 m_bc_values
Definition SSAFEM.hh:68
array::Array2D< Coefficients > m_coefficients
Definition SSAFEM.hh:74
double get_min_thickness() const
Returns minimum thickness to trigger use of extension.
Definition SSA.cc:66
array::Vector m_velocity_global
Definition SSA.hh:122
SSAStrengthExtension * strength_extension
Definition SSA.hh:104
std::shared_ptr< rheology::FlowLaw > m_flow_law
#define PISM_CHK(errcode, name)
#define PISM_ERROR_LOCATION
const int n_chi
Definition FEM.hh:191
const unsigned int MAX_QUADRATURE_SIZE
Definition FEM.hh:173
static double secondInvariant_2D(const Vector2d &U_x, const Vector2d &U_y)
Definition FlowLaw.hh:45
static const double g
Definition exactTestP.cc:36
static const double k
Definition exactTestP.cc:42
double dy
Function derivative with respect to y.
Definition FEM.hh:157
double val
Function value.
Definition FEM.hh:153
double dx
Function deriviative with respect to x.
Definition FEM.hh:155
Struct for gathering the value and derivative of a function at a point.
Definition FEM.hh:151