Loading [MathJax]/extensions/tex2jax.js
PISM, A Parallel Ice Sheet Model 2.2.2-d6b3a29ca committed by Constantine Khrulev on 2025-03-28
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
IPTotalVariationFunctional.cc
Go to the documentation of this file.
1// Copyright (C) 2012, 2013, 2014, 2015, 2016, 2017, 2020, 2022, 2023, 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/functional/IPTotalVariationFunctional.hh"
20#include "pism/util/Grid.hh"
21#include "pism/util/pism_utilities.hh"
22#include "pism/util/array/Scalar.hh"
23#include "pism/util/fem/DirichletData.hh"
24
25namespace pism {
26namespace inverse {
27
29 double c, double exponent, double eps,
30 array::Scalar *dirichletLocations) :
31 IPFunctional<array::Scalar>(grid), m_dirichletIndices(dirichletLocations),
32 m_c(c), m_lebesgue_exp(exponent), m_epsilon_sq(eps*eps) {
33}
34
36
37 const unsigned int Nk = fem::q1::n_chi;
38 const unsigned int Nq = m_element.n_pts();
39 const unsigned int Nq_max = fem::MAX_QUADRATURE_SIZE;
40
41 // The value of the objective
42 double value = 0;
43
44 double x_e[Nk];
45 double x_q[Nq_max], dxdx_q[Nq_max], dxdy_q[Nq_max];
46
47 array::AccessScope list(x);
48
50
51 // Loop through all LOCAL elements.
52 const int
57
58 for (int j = ys; j < ys + ym; j++) {
59 for (int i = xs; i < xs + xm; i++) {
60 m_element.reset(i, j);
61
62 // Obtain values of x at the quadrature points for the element.
64 if (dirichletBC) {
65 dirichletBC.enforce_homogeneous(m_element, x_e);
66 }
67 m_element.evaluate(x_e, x_q, dxdx_q, dxdy_q);
68
69 for (unsigned int q = 0; q < Nq; q++) {
70 auto W = m_element.weight(q);
71 value += m_c*W*pow(m_epsilon_sq + dxdx_q[q]*dxdx_q[q] + dxdy_q[q]*dxdy_q[q], m_lebesgue_exp / 2);
72 } // q
73 } // j
74 } // i
75
76 GlobalSum(m_grid->com, &value, OUTPUT, 1);
77}
78
80
81 const unsigned int Nk = fem::q1::n_chi;
82 const unsigned int Nq = m_element.n_pts();
83 const unsigned int Nq_max = fem::MAX_QUADRATURE_SIZE;
84
85 // Clear the gradient before doing anything with it.
86 gradient.set(0);
87
88 double x_e[Nk];
89 double x_q[Nq_max], dxdx_q[Nq_max], dxdy_q[Nq_max];
90
91 double gradient_e[Nk];
92
93 array::AccessScope list{&x, &gradient};
94
96
97 // Loop through all local and ghosted elements.
98 const int
100 xm = m_element_index.xm,
101 ys = m_element_index.ys,
102 ym = m_element_index.ym;
103
104 for (int j = ys; j < ys + ym; j++) {
105 for (int i = xs; i < xs + xm; i++) {
106
107 // Reset the DOF map for this element.
108 m_element.reset(i, j);
109
110 // Obtain values of x at the quadrature points for the element.
111 m_element.nodal_values(x.array(), x_e);
112 if (dirichletBC) {
113 dirichletBC.constrain(m_element);
114 dirichletBC.enforce_homogeneous(m_element, x_e);
115 }
116 m_element.evaluate(x_e, x_q, dxdx_q, dxdy_q);
117
118 // Zero out the element-local residual in preparation for updating it.
119 for (unsigned int k = 0; k < Nk; k++) {
120 gradient_e[k] = 0;
121 }
122
123 for (unsigned int q = 0; q < Nq; q++) {
124 auto W = m_element.weight(q);
125 const double &dxdx_qq = dxdx_q[q], &dxdy_qq = dxdy_q[q];
126 for (unsigned int k = 0; k < Nk; k++) {
127 gradient_e[k] += m_c*W*(m_lebesgue_exp)*pow(m_epsilon_sq + dxdx_q[q]*dxdx_q[q] + dxdy_q[q]*dxdy_q[q], m_lebesgue_exp / 2 - 1)
128 *(dxdx_qq*m_element.chi(q, k).dx + dxdy_qq*m_element.chi(q, k).dy);
129 } // k
130 } // q
131 m_element.add_contribution(gradient_e, gradient.array());
132 } // j
133 } // i
134}
135
136} // end of namespace inverse
137} // end of namespace pism
Makes sure that we call begin_access() and end_access() for all accessed array::Arrays.
Definition Array.hh:64
void set(double c)
Result: v[j] <- c for all j.
Definition Array.cc:629
void enforce_homogeneous(const Element2 &element, double *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 lym
total number local elements in y direction.
int lxm
total number local elements in x direction.
int lxs
x-index of the first local element.
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 lys
y-index of the first local element.
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
std::shared_ptr< const Grid > m_grid
Abstract base class for functions from ice model vectors to .
virtual void valueAt(array::Scalar &x, double *OUTPUT)
Computes the value of the functional at the vector x.
virtual void gradientAt(array::Scalar &x, array::Scalar &gradient)
Computes the gradient of the functional at the vector x.
IPTotalVariationFunctional2S(std::shared_ptr< const Grid > grid, double c, double q, double eps, array::Scalar *dirichletLocations=NULL)
const int n_chi
Definition FEM.hh:191
const unsigned int MAX_QUADRATURE_SIZE
Definition FEM.hh:173
static const double k
Definition exactTestP.cc:42
void GlobalSum(MPI_Comm comm, double *local, double *result, int count)
double dy
Function derivative with respect to y.
Definition FEM.hh:157
double dx
Function deriviative with respect to x.
Definition FEM.hh:155