PISM, A Parallel Ice Sheet Model  stable v2.1-1-g6902d5502 committed by Ed Bueler on 2023-12-20 08:38:27 -0800
IP_SSATaucForwardProblem.cc
Go to the documentation of this file.
1 // Copyright (C) 2012, 2014, 2015, 2016, 2017, 2019, 2020, 2021, 2022, 2023 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_SSATaucForwardProblem.hh"
20 #include "pism/basalstrength/basal_resistance.hh"
21 #include "pism/util/Grid.hh"
22 #include "pism/util/Mask.hh"
23 #include "pism/util/Vars.hh"
24 #include "pism/util/error_handling.hh"
25 #include "pism/util/pism_utilities.hh"
26 #include "pism/geometry/Geometry.hh"
27 #include "pism/stressbalance/StressBalance.hh"
28 #include "pism/util/petscwrappers/DM.hh"
29 #include "pism/util/petscwrappers/Vec.hh"
30 
31 namespace pism {
32 namespace inverse {
33 
34 IP_SSATaucForwardProblem::IP_SSATaucForwardProblem(std::shared_ptr<const Grid> grid,
36  : SSAFEM(grid),
37  m_zeta(NULL),
38  m_dzeta_local(m_grid, "d_zeta_local"),
39  m_tauc_copy(m_grid, "tauc"),
40  m_fixed_tauc_locations(NULL),
41  m_tauc_param(tp),
42  m_du_global(m_grid, "linearization work vector (sans ghosts)"),
43  m_du_local(m_grid, "linearization work vector (with ghosts)"),
44  m_element_index(*m_grid),
45  m_element(*m_grid, fem::Q1Quadrature4()),
46  m_rebuild_J_state(true)
47 {
48 
49  PetscErrorCode ierr;
50 
51  m_velocity_shared.reset(new array::Vector(m_grid, "dummy"));
52  m_velocity_shared->metadata(0) = m_velocity.metadata(0);
53  m_velocity_shared->metadata(1) = m_velocity.metadata(1);
54 
56  .long_name("yield stress for basal till (plastic or pseudo-plastic model)")
57  .units("Pa");
58 
59  ierr = DMSetMatType(*m_da, MATBAIJ);
60  PISM_CHK(ierr, "DMSetMatType");
61  ierr = DMCreateMatrix(*m_da, m_J_state.rawptr());
62  PISM_CHK(ierr, "DMCreateMatrix");
63 
64  ierr = KSPCreate(m_grid->com, m_ksp.rawptr());
65  PISM_CHK(ierr, "KSPCreate");
66 
67  double ksp_rtol = 1e-12;
68  ierr = KSPSetTolerances(m_ksp, ksp_rtol, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT);
69  PISM_CHK(ierr, "KSPSetTolerances");
70 
71  PC pc;
72  ierr = KSPGetPC(m_ksp, &pc);
73  PISM_CHK(ierr, "KSPGetPC");
74 
75  ierr = PCSetType(pc, PCBJACOBI);
76  PISM_CHK(ierr, "PCSetType");
77 
78  ierr = KSPSetFromOptions(m_ksp);
79  PISM_CHK(ierr, "KSPSetFromOptions");
80 }
81 
83 
84  // This calls SSA::init(), which calls pism::Vars::get_2d_scalar()
85  // to set m_tauc.
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);
96  geometry.ice_area_specific_volume.set(0.0);
97 
98  geometry.ensure_consistency(m_config->get_number("stress_balance.ice_free_thickness_standard"));
99 
100  stressbalance::Inputs inputs;
101 
102  inputs.geometry = &geometry;
103  inputs.basal_melt_rate = NULL;
104  inputs.basal_yield_stress = m_grid->variables().get_2d_scalar("tauc");
105  inputs.enthalpy = m_grid->variables().get_3d_scalar("enthalpy");
106  inputs.age = NULL;
107  inputs.bc_mask = m_grid->variables().get_2d_mask("vel_bc_mask");
108  inputs.bc_values = m_grid->variables().get_2d_vector("vel_bc");
109 
110  inputs.water_column_pressure = NULL;
111 
112  cache_inputs(inputs);
113  }
114 }
115 
116 //! Sets the current value of of the design paramter \f$\zeta\f$.
117 /*! This method sets \f$\zeta\f$ but does not solve the %SSA.
118 It it intended for inverse methods that simultaneously compute
119 the pair \f$u\f$ and \f$\zeta\f$ without ever solving the %SSA
120 directly. Use this method in conjuction with
121 \ref assemble_jacobian_state and \ref apply_jacobian_design and their friends.
122 The vector \f$\zeta\f$ is not copied; a reference to the array::Array is
123 kept.
124 */
126 
127  array::Scalar &tauc = m_tauc_copy;
128 
129  m_zeta = &new_zeta;
130 
131  // Convert zeta to tauc.
133 
134  // Cache tauc at the quadrature points.
135  array::AccessScope list{&tauc, &m_coefficients};
136 
137  for (auto p = m_grid->points(1); p; p.next()) {
138  const int i = p.i(), j = p.j();
139  m_coefficients(i, j).tauc = tauc(i, j);
140  }
141 
142  // Flag the state jacobian as needing rebuilding.
143  m_rebuild_J_state = true;
144 }
145 
146 //! Sets the current value of the design variable \f$\zeta\f$ and solves the %SSA to find the associated \f$u_{\rm SSA}\f$.
147 /* Use this method for inverse methods employing the reduced gradient. Use this method
148 in conjuction with apply_linearization and apply_linearization_transpose.*/
149 std::shared_ptr<TerminationReason> IP_SSATaucForwardProblem::linearize_at(array::Scalar &zeta) {
150  this->set_design(zeta);
151  return this->solve_nocache();
152 }
153 
154 //! Computes the residual function \f$\mathcal{R}(u, \zeta)\f$ as defined in the class-level documentation.
155 /* The value of \f$\zeta\f$ is set prior to this call via set_design or linearize_at. The value
156 of the residual is returned in \a RHS.*/
158  array::AccessScope l{&u, &RHS};
159 
160  this->compute_local_function(u.array(), RHS.array());
161 }
162 
163 //! Computes the residual function \f$\mathcal{R}(u, \zeta)\f$ defined in the class-level documentation.
164 /* The return value is specified via a Vec for the benefit of certain TAO routines. Otherwise,
165 the method is identical to the assemble_residual returning values as a StateVec (an array::Vector).*/
167  array::AccessScope l{&u};
168 
169  petsc::DMDAVecArray rhs_a(m_da, RHS);
170  this->compute_local_function(u.array(), (Vector2d**)rhs_a.get());
171 }
172 
173 //! Assembles the state Jacobian matrix.
174 /* The matrix depends on the current value of the design variable \f$\zeta\f$ and the current
175 value of the state variable \f$u\f$. The specification of \f$\zeta\f$ is done earlier
176 with set_design or linearize_at. The value of \f$u\f$ is specified explicitly as an argument
177 to this method.
178  @param[in] u Current state variable value.
179  @param[out] J computed state Jacobian.
180 */
182  array::AccessScope l{&u};
183 
184  this->compute_local_jacobian(u.array(), Jac);
185 }
186 
187 //! Applies the design Jacobian matrix to a perturbation of the design variable.
188 /*! The return value uses a DesignVector (array::Vector), which can be ghostless. Ghosts (if present) are updated.
189 \overload
190 */
192  array::Vector &du) {
193  array::AccessScope l{&du};
194 
195  this->apply_jacobian_design(u, dzeta, du.array());
196 }
197 
198 //! Applies the design Jacobian matrix to a perturbation of the design variable.
199 /*! The return value is a Vec for the benefit of TAO. It is assumed to
200 be ghostless; no communication is done. \overload
201 */
203  Vec du) {
204  petsc::DMDAVecArray du_a(m_da, du);
205  this->apply_jacobian_design(u, dzeta, (Vector2d**)du_a.get());
206 }
207 
208 //! Applies the design Jacobian matrix to a perturbation of the design variable.
209 /*! The matrix depends on the current value of the design variable \f$\zeta\f$ and the current
210 value of the state variable \f$u\f$. The specification of \f$\zeta\f$ is done earlier
211 with set_design or linearize_at. The value of \f$u\f$ is specified explicitly as an argument
212 to this method.
213  @param[in] u Current state variable value.
214  @param[in] dzeta Perturbation of the design variable. Prefers vectors with ghosts; will copy to a ghosted vector if needed.
215  @param[out] du_a Computed corresponding perturbation of the state variable. The array \a du_a
216  should be extracted first from a Vec or an array::Array.
217 
218  Typically this method is called via one of its overloads.
219 */
221  array::Scalar &dzeta,
222  Vector2d **du_a) {
223  const unsigned int Nk = fem::q1::n_chi;
224  const unsigned int Nq = m_element.n_pts();
225  const unsigned int Nq_max = fem::MAX_QUADRATURE_SIZE;
226 
228 
229  array::Scalar *dzeta_local;
230  if (dzeta.stencil_width() > 0) {
231  dzeta_local = &dzeta;
232  } else {
233  m_dzeta_local.copy_from(dzeta);
234  dzeta_local = &m_dzeta_local;
235  }
236  list.add(*dzeta_local);
237 
238  // Zero out the portion of the function we are responsible for computing.
239  for (auto p = m_grid->points(); p; p.next()) {
240  const int i = p.i(), j = p.j();
241 
242  du_a[j][i].u = 0.0;
243  du_a[j][i].v = 0.0;
244  }
245 
246  // Aliases to help with notation consistency below.
247  const array::Scalar *dirichletLocations = &m_bc_mask;
248  const array::Vector *dirichletValues = &m_bc_values;
249  double dirichletWeight = m_dirichletScale;
250 
251  Vector2d u_e[Nk];
252  Vector2d u_q[Nq_max];
253 
254  Vector2d du_e[Nk];
255 
256  double dzeta_e[Nk];
257 
258  double zeta_e[Nk];
259 
260  double dtauc_e[Nk];
261  double dtauc_q[Nq_max];
262 
263  fem::DirichletData_Vector dirichletBC(dirichletLocations, dirichletValues,
264  dirichletWeight);
266 
267  // Loop through all elements.
268  const int
269  xs = m_element_index.xs,
270  xm = m_element_index.xm,
271  ys = m_element_index.ys,
272  ym = m_element_index.ym;
273 
274  ParallelSection loop(m_grid->com);
275  try {
276  for (int j = ys; j < ys + ym; j++) {
277  for (int i = xs; i < xs + xm; i++) {
278 
279  // Zero out the element residual in prep for updating it.
280  for (unsigned int k = 0; k < Nk; k++) {
281  du_e[k].u = 0;
282  du_e[k].v = 0;
283  }
284 
285  // Initialize the map from global to local degrees of freedom for this element.
286  m_element.reset(i, j);
287 
288  // Obtain the value of the solution at the nodes adjacent to the element,
289  // fix dirichlet values, and compute values at quad pts.
290  m_element.nodal_values(u.array(), u_e);
291  if (dirichletBC) {
292  dirichletBC.constrain(m_element);
293  dirichletBC.enforce(m_element, u_e);
294  }
295  m_element.evaluate(u_e, u_q);
296 
297  // Compute dzeta at the nodes
298  m_element.nodal_values(dzeta_local->array(), dzeta_e);
299  if (fixedZeta) {
300  fixedZeta.enforce_homogeneous(m_element, dzeta_e);
301  }
302 
303  // Compute the change in tau_c with respect to zeta at the quad points.
304  m_element.nodal_values(m_zeta->array(), zeta_e);
305  for (unsigned int k = 0; k < Nk; k++) {
306  m_tauc_param.toDesignVariable(zeta_e[k], NULL, dtauc_e + k);
307  dtauc_e[k] *= dzeta_e[k];
308  }
309  m_element.evaluate(dtauc_e, dtauc_q);
310 
311  int mask[Nq_max];
312  {
313  Coefficients coeffs[Nk];
314  double thickness[Nq_max];
315  double tauc[Nq_max];
316  double hardness[Nq_max];
317 
318  m_element.nodal_values(m_coefficients.array(), coeffs);
319 
321  mask, thickness, tauc, hardness);
322  }
323 
324  for (unsigned int q = 0; q < Nq; q++) {
325  Vector2d u_qq = u_q[q];
326 
327  // Determine "dbeta / dzeta" at the quadrature point
328  double dbeta = 0;
329  if (mask::grounded_ice(mask[q])) {
330  dbeta = m_basal_sliding_law->drag(dtauc_q[q], u_qq.u, u_qq.v);
331  }
332 
333  auto W = m_element.weight(q);
334 
335  for (unsigned int k = 0; k < Nk; k++) {
336  du_e[k].u += W*dbeta*u_qq.u*m_element.chi(q, k).val;
337  du_e[k].v += W*dbeta*u_qq.v*m_element.chi(q, k).val;
338  }
339  } // q
340  m_element.add_contribution(du_e, du_a);
341  } // j
342  } // i
343  } catch (...) {
344  loop.failed();
345  }
346  loop.check();
347 
348  if (dirichletBC) {
349  dirichletBC.fix_residual_homogeneous(du_a);
350  }
351 }
352 
353 //! Applies the transpose of the design Jacobian matrix to a perturbation of the state variable.
354 /*! The return value uses a StateVector (array::Scalar) which can be ghostless; ghosts (if present) are updated.
355 \overload
356 */
358  array::Vector &du,
359  array::Scalar &dzeta) {
360  array::AccessScope l{&dzeta};
361 
362  this->apply_jacobian_design_transpose(u, du, dzeta.array());
363 }
364 
365 //! Applies the transpose of the design Jacobian matrix to a perturbation of the state variable.
366 /*! The return value uses a Vec for the benefit of TAO. It is assumed to be ghostless; no communication is done.
367 \overload */
369  array::Vector &du,
370  Vec dzeta) {
371  petsc::DM::Ptr da2 = m_grid->get_dm(1, m_config->get_number("grid.max_stencil_width"));
372  petsc::DMDAVecArray dzeta_a(da2, dzeta);
373  this->apply_jacobian_design_transpose(u, du, (double**)dzeta_a.get());
374 }
375 
376 //! Applies the transpose of the design Jacobian matrix to a perturbation of the state variable.
377 /*! The matrix depends on the current value of the design variable \f$\zeta\f$ and the current
378 value of the state variable \f$u\f$. The specification of \f$\zeta\f$ is done earlier
379 with set_design or linearize_at. The value of \f$u\f$ is specified explicitly as an argument
380 to this method.
381  @param[in] u Current state variable value.
382  @param[in] du Perturbation of the state variable. Prefers vectors with ghosts; will copy to a ghosted vector if need be.
383  @param[out] dzeta_a Computed corresponding perturbation of the design variable. The array \a dzeta_a
384  should be extracted first from a Vec or an array::Array.
385 
386  Typically this method is called via one of its overloads.
387 */
389  array::Vector &du,
390  double **dzeta_a) {
391  const unsigned int Nk = fem::q1::n_chi;
392  const unsigned int Nq = m_element.n_pts();
393  const unsigned int Nq_max = fem::MAX_QUADRATURE_SIZE;
394 
396 
397  array::Vector *du_local;
398  if (du.stencil_width() > 0) {
399  du_local = &du;
400  } else {
401  m_du_local.copy_from(du);
402  du_local = &m_du_local;
403  }
404  list.add(*du_local);
405 
406  Vector2d u_e[Nk];
407  Vector2d u_q[Nq_max];
408 
409  Vector2d du_e[Nk];
410  Vector2d du_q[Nq_max];
411 
412  double dzeta_e[Nk];
413 
414  // Aliases to help with notation consistency.
415  const array::Scalar *dirichletLocations = &m_bc_mask;
416  const array::Vector *dirichletValues = &m_bc_values;
417  double dirichletWeight = m_dirichletScale;
418 
419  fem::DirichletData_Vector dirichletBC(dirichletLocations, dirichletValues,
420  dirichletWeight);
421 
422  // Zero out the portion of the function we are responsible for computing.
423  for (auto p = m_grid->points(); p; p.next()) {
424  const int i = p.i(), j = p.j();
425 
426  dzeta_a[j][i] = 0;
427  }
428 
429  const int
430  xs = m_element_index.xs,
431  xm = m_element_index.xm,
432  ys = m_element_index.ys,
433  ym = m_element_index.ym;
434 
435  ParallelSection loop(m_grid->com);
436  try {
437  for (int j=ys; j<ys+ym; j++) {
438  for (int i=xs; i<xs+xm; i++) {
439  // Initialize the map from global to local degrees of freedom for this element.
440  m_element.reset(i, j);
441 
442  // Obtain the value of the solution at the nodes adjacent to the element.
443  // Compute the solution values and symmetric gradient at the quadrature points.
444  m_element.nodal_values(du_local->array(), du_e);
445  if (dirichletBC) {
446  dirichletBC.enforce_homogeneous(m_element, du_e);
447  }
448  m_element.evaluate(du_e, du_q);
449 
450  m_element.nodal_values(u.array(), u_e);
451  if (dirichletBC) {
452  dirichletBC.enforce(m_element, u_e);
453  }
454  m_element.evaluate(u_e, u_q);
455 
456  // Zero out the element-local residual in prep for updating it.
457  for (unsigned int k=0; k<Nk; k++) {
458  dzeta_e[k] = 0;
459  }
460 
461  int mask[Nq_max];
462  {
463  Coefficients coeffs[Nk];
464  double thickness[Nq_max];
465  double tauc[Nq_max];
466  double hardness[Nq_max];
467 
468  m_element.nodal_values(m_coefficients.array(), coeffs);
469 
471  mask, thickness, tauc, hardness);
472  }
473 
474  for (unsigned int q=0; q<Nq; q++) {
475  Vector2d du_qq = du_q[q];
476  Vector2d u_qq = u_q[q];
477 
478  // Determine "dbeta/dtauc" at the quadrature point
479  double dbeta_dtauc = 0;
480  if (mask::grounded_ice(mask[q])) {
481  dbeta_dtauc = m_basal_sliding_law->drag(1., u_qq.u, u_qq.v);
482  }
483 
484  auto W = m_element.weight(q);
485 
486  for (unsigned int k=0; k<Nk; k++) {
487  dzeta_e[k] += W*dbeta_dtauc*(du_qq.u*u_qq.u+du_qq.v*u_qq.v)*m_element.chi(q, k).val;
488  }
489  } // q
490 
491  m_element.add_contribution(dzeta_e, dzeta_a);
492  } // j
493  } // i
494  } catch (...) {
495  loop.failed();
496  }
497  loop.check();
498 
499  for (auto p = m_grid->points(); p; p.next()) {
500  const int i = p.i(), j = p.j();
501 
502  double dtauc_dzeta;
503  m_tauc_param.toDesignVariable((*m_zeta)(i, j), NULL, &dtauc_dzeta);
504  dzeta_a[j][i] *= dtauc_dzeta;
505  }
506 
509  fixedTauc.fix_residual_homogeneous(dzeta_a);
510  }
511 }
512 
513 /*!\brief Applies the linearization of the forward map (i.e. the reduced gradient \f$DF\f$ described in
514 the class-level documentation.) */
515 /*! As described previously,
516 \f[
517 Df = J_{\rm State}^{-1} J_{\rm Design}.
518 \f]
519 Applying the linearization then involves the solution of a linear equation.
520 The matrices \f$J_{\rm State}\f$ and \f$J_{\rm Design}\f$ both depend on the value of the
521 design variable \f$\zeta\f$ and the value of the corresponding state variable \f$u=F(\zeta)\f$.
522 These are established by first calling linearize_at.
523  @param[in] dzeta Perturbation of the design variable
524  @param[out] du Computed corresponding perturbation of the state variable; ghosts (if present) are updated.
525 */
527 
528  PetscErrorCode ierr;
529 
530  if (m_rebuild_J_state) {
532  m_rebuild_J_state = false;
533  }
534 
536  m_du_global.scale(-1);
537 
538  // call PETSc to solve linear system by iterative method.
539  ierr = KSPSetOperators(m_ksp, m_J_state, m_J_state);
540  PISM_CHK(ierr, "KSPSetOperators");
541 
542  ierr = KSPSolve(m_ksp, m_du_global.vec(), m_du_global.vec());
543  PISM_CHK(ierr, "KSPSolve"); // SOLVE
544 
545  KSPConvergedReason reason;
546  ierr = KSPGetConvergedReason(m_ksp, &reason);
547  PISM_CHK(ierr, "KSPGetConvergedReason");
548  if (reason < 0) {
550  "IP_SSATaucForwardProblem::apply_linearization solve"
551  " failed to converge (KSP reason %s)",
552  KSPConvergedReasons[reason]);
553  } else {
554  m_log->message(4,
555  "IP_SSATaucForwardProblem::apply_linearization converged"
556  " (KSP reason %s)\n",
557  KSPConvergedReasons[reason]);
558  }
559 
561 }
562 
563 /*! \brief Applies the transpose of the linearization of the forward map
564  (i.e. the transpose of the reduced gradient \f$DF\f$ described in the class-level documentation.) */
565 /*! As described previously,
566 \f[
567 Df = J_{\rm State}^{-1} J_{\rm Design}.
568 \f]
569 so
570 \f[
571 Df^t = J_{\rm Design}^t \; (J_{\rm State}^t)^{-1} .
572 \f]
573 Applying the transpose of the linearization then involves the solution of a linear equation.
574 The matrices \f$J_{\rm State}\f$ and \f$J_{\rm Design}\f$ both depend on the value of the
575 design variable \f$\zeta\f$ and the value of the corresponding state variable \f$u=F(\zeta)\f$.
576 These are established by first calling linearize_at.
577  @param[in] du Perturbation of the state variable
578  @param[out] dzeta Computed corresponding perturbation of the design variable; ghosts (if present) are updated.
579 */
581  array::Scalar &dzeta) {
582 
583  PetscErrorCode ierr;
584 
585  if (m_rebuild_J_state) {
587  m_rebuild_J_state = false;
588  }
589 
590  // Aliases to help with notation consistency below.
591  const array::Scalar *dirichletLocations = &m_bc_mask;
592  const array::Vector *dirichletValues = &m_bc_values;
593  double dirichletWeight = m_dirichletScale;
594 
596 
597  fem::DirichletData_Vector dirichletBC(dirichletLocations, dirichletValues, dirichletWeight);
598 
599  if (dirichletBC) {
602  }
603 
604  // call PETSc to solve linear system by iterative method.
605  ierr = KSPSetOperators(m_ksp, m_J_state, m_J_state);
606  PISM_CHK(ierr, "KSPSetOperators");
607 
608  ierr = KSPSolve(m_ksp, m_du_global.vec(), m_du_global.vec());
609  PISM_CHK(ierr, "KSPSolve"); // SOLVE
610 
611  KSPConvergedReason reason;
612  ierr = KSPGetConvergedReason(m_ksp, &reason);
613  PISM_CHK(ierr, "KSPGetConvergedReason");
614 
615  if (reason < 0) {
617  "IP_SSATaucForwardProblem::apply_linearization solve"
618  " failed to converge (KSP reason %s)",
619  KSPConvergedReasons[reason]);
620  } else {
621  m_log->message(4,
622  "IP_SSATaucForwardProblem::apply_linearization converged"
623  " (KSP reason %s)\n",
624  KSPConvergedReasons[reason]);
625  }
626 
628  dzeta.scale(-1);
629 
630  if (dzeta.stencil_width() > 0) {
631  dzeta.update_ghosts();
632  }
633 }
634 
635 } // end of namespace inverse
636 } // end of namespace pism
const Config::ConstPtr m_config
configuration database used by this component
Definition: Component.hh:158
const Logger::ConstPtr m_log
logger (for easy access)
Definition: Component.hh:162
const std::shared_ptr< const Grid > m_grid
grid used by this component
Definition: Component.hh:156
array::Scalar1 sea_level_elevation
Definition: Geometry.hh:48
void ensure_consistency(double ice_free_thickness_threshold)
Definition: Geometry.cc:112
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
virtual double drag(double tauc, double vx, double vy) const
Compute the drag coefficient for the basal shear stress.
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
VariableMetadata & long_name(const std::string &input)
VariableMetadata & units(const std::string &input)
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:39
std::shared_ptr< Wrapper > Ptr
Definition: Wrapper.hh:30
Makes sure that we call begin_access() and end_access() for all accessed array::Arrays.
Definition: Array.hh:65
void copy_from(const Array2D< T > &source)
Definition: Array2D.hh:73
void add(double alpha, const Array2D< T > &x)
Definition: Array2D.hh:65
petsc::Vec & vec() const
Definition: Array.cc:339
void scale(double alpha)
Result: v <- v * alpha. Calls VecScale.
Definition: Array.cc:253
void set(double c)
Result: v[j] <- c for all j.
Definition: Array.cc:707
void update_ghosts()
Updates ghost points.
Definition: Array.cc:693
unsigned int stencil_width() const
Get the stencil width of the current Array. Returns 0 if ghosts are not available.
Definition: Array.cc:331
SpatialVariableMetadata & metadata(unsigned int N=0)
Returns a reference to the SpatialVariableMetadata object containing metadata for the compoment N.
Definition: Array.cc:553
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:187
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:231
void nodal_values(const array::Scalar &x_global, int *result) const
Get nodal values of an integer mask.
Definition: Element.cc:176
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.
const Germ & chi(unsigned int q, unsigned int k) const
Definition: Element.hh:73
int n_pts() const
Number of quadrature points.
Definition: Element.hh:80
double weight(unsigned int q) const
Weight of the quadrature point q
Definition: Element.hh:85
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 .
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 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::Scalar * m_fixed_tauc_locations
Locations where should not be adjusted.
IPDesignVariableParameterization & m_tauc_param
The function taking to .
bool m_rebuild_J_state
Flag indicating that the state jacobian matrix needs rebuilding.
array::Scalar * m_zeta
Current value of zeta, provided from caller.
petsc::KSP m_ksp
KSP used in apply_linearization and apply_linearization_transpose.
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...
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 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 .
std::shared_ptr< array::Vector > m_velocity_shared
Copy of the velocity field managed using a shared pointer.
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...
array::Vector1 m_du_local
Temporary storage when state vectors need to be used with ghosts.
array::Vector m_du_global
Temporary storage when state vectors need to be used without ghosts.
IP_SSATaucForwardProblem(std::shared_ptr< const Grid > g, IPDesignVariableParameterization &tp)
virtual void set_design(array::Scalar &zeta)
Sets the current value of of the design paramter .
array::Scalar1 m_dzeta_local
Storage for d_zeta with ghosts, if needed when an argument d_zeta is ghost-less.
virtual void assemble_jacobian_state(array::Vector &u, Mat J)
Assembles the state Jacobian matrix.
array::Scalar2 m_tauc_copy
Storage for tauc (avoids modifying fields obtained via pism::Vars)
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 compute_local_function(Vector2d const *const *const velocity, Vector2d **residual)
Implements the callback for computing the residual.
Definition: SSAFEM.cc:726
void quad_point_values(const fem::Element &Q, 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:336
void cache_inputs(const Inputs &inputs)
Initialize stored data from the coefficients in the SSA. Called by SSAFEM::solve.
Definition: SSAFEM.cc:258
void compute_local_jacobian(Vector2d const *const *const velocity, Mat J)
Implements the callback for computing the Jacobian.
Definition: SSAFEM.cc:948
array::Scalar1 m_bc_mask
Definition: SSAFEM.hh:69
std::shared_ptr< TerminationReason > solve_nocache()
Definition: SSAFEM.cc:184
array::Vector1 m_bc_values
Definition: SSAFEM.hh:70
array::Array2D< Coefficients > m_coefficients
Definition: SSAFEM.hh:76
std::shared_ptr< petsc::DM > m_da
Definition: SSA.hh:149
IceBasalResistancePlasticLaw * m_basal_sliding_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
bool grounded_ice(int M)
Definition: Mask.hh:51
static const double k
Definition: exactTestP.cc:42
double val
Function value.
Definition: FEM.hh:153