PISM, A Parallel Ice Sheet Model 2.3.0-79cae578d committed by Constantine Khrulev on 2026-03-22
Loading...
Searching...
No Matches
IPGroundedIceH1NormFunctional.cc
Go to the documentation of this file.
1// Copyright (C) 2013, 2014, 2015, 2016, 2017, 2020, 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/IPGroundedIceH1NormFunctional.hh"
20#include "pism/util/error_handling.hh"
21#include "pism/util/Grid.hh"
22#include "pism/util/pism_utilities.hh"
23#include "pism/util/array/CellType.hh"
24#include "pism/util/fem/DirichletData.hh"
25
26namespace pism {
27namespace inverse {
28
30
31 const unsigned int Nk = fem::q1::n_chi;
32 const unsigned int Nq = m_element.n_pts();
33 const unsigned int Nq_max = fem::MAX_QUADRATURE_SIZE;
34
35 // The value of the objective
36 double value = 0;
37
38 double x_e[Nk];
39 double x_q[Nq_max], dxdx_q[Nq_max], dxdy_q[Nq_max];
40
41 // initialize x_e to silence the static analyzer
42 for (unsigned int k = 0; k < Nk; ++k) {
43 x_e[k] = 0.0;
44 }
45
47
49
50 // Loop through all LOCAL elements.
51 const int
56
57 for (int j=ys; j<ys+ym; j++) {
58 for (int i=xs; i<xs+xm; i++) {
59 bool all_grounded_ice = (m_ice_mask.grounded_ice(i, j) and
60 m_ice_mask.grounded_ice(i+1, j) and
61 m_ice_mask.grounded_ice(i, j+1) and
62 m_ice_mask.grounded_ice(i+1, j+1));
63
64 if (! all_grounded_ice) {
65 continue;
66 }
67
68 m_element.reset(i, j);
69
70 // Obtain values of x at the quadrature points for the element.
72 if (dirichletBC) {
73 dirichletBC.enforce_homogeneous(m_element, x_e);
74 }
75 m_element.evaluate(x_e, x_q, dxdx_q, dxdy_q);
76
77 for (unsigned int q=0; q<Nq; q++) {
78 auto W = m_element.weight(q);
79 value += W*(m_cL2*x_q[q]*x_q[q]+ m_cH1*(dxdx_q[q]*dxdx_q[q]+dxdy_q[q]*dxdy_q[q]));
80 } // q
81 } // j
82 } // i
83
84 GlobalSum(m_grid->com, &value, OUTPUT, 1);
85}
86
88
89 const unsigned int Nk = fem::q1::n_chi;
90 const unsigned int Nq = m_element.n_pts();
91 const unsigned int Nq_max = fem::MAX_QUADRATURE_SIZE;
92
93 // The value of the objective
94 double value = 0;
95
96 double a_e[Nk];
97 double a_q[Nq_max], dadx_q[Nq_max], dady_q[Nq_max];
98
99 array::AccessScope list{&a, &b, &m_ice_mask};
100
101 double b_e[Nk];
102 double b_q[Nq_max], dbdx_q[Nq_max], dbdy_q[Nq_max];
103
104 // initialize a_e and b_e to silence the static analyzer
105 for (unsigned int k = 0; k < Nk; ++k) {
106 a_e[k] = 0.0;
107 b_e[k] = 0.0;
108 }
109
111
112 // Loop through all LOCAL elements.
113 const int
114 xs = m_element_index.lxs,
115 xm = m_element_index.lxm,
116 ys = m_element_index.lys,
117 ym = m_element_index.lym;
118
119 for (int j=ys; j<ys+ym; j++) {
120 for (int i=xs; i<xs+xm; i++) {
121 bool all_grounded_ice = (m_ice_mask.grounded_ice(i, j) and
122 m_ice_mask.grounded_ice(i+1, j) and
123 m_ice_mask.grounded_ice(i, j+1) and
124 m_ice_mask.grounded_ice(i+1, j+1));
125
126 if (! all_grounded_ice) {
127 continue;
128 }
129
130 m_element.reset(i, j);
131
132 // Obtain values of x at the quadrature points for the element.
133 m_element.nodal_values(a.array(), a_e);
134 if (dirichletBC) {
135 dirichletBC.enforce_homogeneous(m_element, a_e);
136 }
137 m_element.evaluate(a_e, a_q, dadx_q, dady_q);
138
139 m_element.nodal_values(b.array(), b_e);
140 if (dirichletBC) {
141 dirichletBC.enforce_homogeneous(m_element, b_e);
142 }
143 m_element.evaluate(b_e, b_q, dbdx_q, dbdy_q);
144
145 for (unsigned int q=0; q<Nq; q++) {
146 auto W = m_element.weight(q);
147 value += W*(m_cL2*a_q[q]*b_q[q]+ m_cH1*(dadx_q[q]*dbdx_q[q]+dady_q[q]*dbdy_q[q]));
148 } // q
149 } // j
150 } // i
151
152 GlobalSum(m_grid->com, &value, OUTPUT, 1);
153}
154
155
157
158 const unsigned int Nk = fem::q1::n_chi;
159 const unsigned int Nq = m_element.n_pts();
160 const unsigned int Nq_max = fem::MAX_QUADRATURE_SIZE;
161
162 // Clear the gradient before doing anything with it!
163 gradient.set(0);
164
165 double x_e[Nk];
166 double x_q[Nq_max], dxdx_q[Nq_max], dxdy_q[Nq_max];
167
168 // initialize x_e to silence the static analyzer
169 for (unsigned int k = 0; k < Nk; ++k) {
170 x_e[k] = 0.0;
171 }
172
173 array::AccessScope list{&x, &gradient, &m_ice_mask};
174
175 double gradient_e[Nk];
176
178
179 // Loop through all local and ghosted elements.
180 const int
181 xs = m_element_index.xs,
182 xm = m_element_index.xm,
183 ys = m_element_index.ys,
184 ym = m_element_index.ym;
185
186 for (int j=ys; j<ys+ym; j++) {
187 for (int i=xs; i<xs+xm; i++) {
188 bool all_grounded_ice = (m_ice_mask.grounded_ice(i, j) and
189 m_ice_mask.grounded_ice(i+1, j) and
190 m_ice_mask.grounded_ice(i, j+1) and
191 m_ice_mask.grounded_ice(i+1, j+1));
192
193 if (! all_grounded_ice) {
194 continue;
195 }
196
197 // Reset the DOF map for this element.
198 m_element.reset(i, j);
199
200 // Obtain values of x at the quadrature points for the element.
201 m_element.nodal_values(x.array(), x_e);
202 if (dirichletBC) {
203 dirichletBC.constrain(m_element);
204 dirichletBC.enforce_homogeneous(m_element, x_e);
205 }
206 m_element.evaluate(x_e, x_q, dxdx_q, dxdy_q);
207
208 // Zero out the element-local residual in prep for updating it.
209 for (unsigned int k=0; k<Nk; k++) {
210 gradient_e[k] = 0;
211 }
212
213 for (unsigned int q=0; q<Nq; q++) {
214 auto W = m_element.weight(q);
215 const double &x_qq=x_q[q];
216 const double &dxdx_qq=dxdx_q[q], &dxdy_qq=dxdy_q[q];
217 for (unsigned int k=0; k<Nk; k++) {
218 gradient_e[k] += 2*W*(m_cL2*x_qq*m_element.chi(q, k).val +
219 m_cH1*(dxdx_qq*m_element.chi(q, k).dx + dxdy_qq*m_element.chi(q, k).dy));
220 } // k
221 } // q
222 m_element.add_contribution(gradient_e, gradient.array());
223 } // j
224 } // i
225}
226
228
229 const unsigned int Nk = fem::q1::n_chi;
230 const unsigned int Nq = m_element.n_pts();
231
232 PetscErrorCode ierr;
233 int i, j;
234
235 // Zero out the Jacobian in preparation for updating it.
236 ierr = MatZeroEntries(form);
237 PISM_CHK(ierr, "MatZeroEntries");
238
240
242
243 // Loop through all the elements.
244 const int
245 xs = m_element_index.xs,
246 xm = m_element_index.xm,
247 ys = m_element_index.ys,
248 ym = m_element_index.ym;
249
250 for (j=ys; j<ys+ym; j++) {
251 for (i=xs; i<xs+xm; i++) {
252 bool all_grounded_ice = (m_ice_mask.grounded_ice(i, j) and
253 m_ice_mask.grounded_ice(i+1, j) and
254 m_ice_mask.grounded_ice(i, j+1) and
255 m_ice_mask.grounded_ice(i+1, j+1));
256 if (! all_grounded_ice) {
257 continue;
258 }
259
260 // Element-local Jacobian matrix (there are Nk vector valued degrees
261 // of freedom per elment, for a total of (2*Nk)*(2*Nk) = 16
262 // entries in the local Jacobian.
263 double K[Nk][Nk];
264
265 // Initialize the map from global to local degrees of freedom for this element.
266 m_element.reset(i, j);
267
268 // Don't update rows/cols where we project to zero.
269 if (zeroLocs) {
270 zeroLocs.constrain(m_element);
271 }
272
273 // Build the element-local Jacobian.
274 ierr = PetscMemzero(K, sizeof(K));
275 PISM_CHK(ierr, "PetscMemzero");
276
277 for (unsigned int q=0; q<Nq; q++) {
278 auto W = m_element.weight(q);
279 for (unsigned int k = 0; k < Nk; k++) { // Test functions
280 const fem::Germ &test_qk=m_element.chi(q, k);
281 for (unsigned int l = 0; l < Nk; l++) { // Trial functions
282 const fem::Germ &test_ql=m_element.chi(q, l);
283 K[k][l] += W*(m_cL2*test_qk.val*test_ql.val
284 + m_cH1*(test_qk.dx*test_ql.dx + test_qk.dy*test_ql.dy));
285 } // l
286 } // k
287 } // q
288 m_element.add_contribution(&K[0][0], form);
289 } // j
290 } // i
291
292 if (zeroLocs) {
293 zeroLocs.fix_jacobian(form);
294 }
295
296 ierr = MatAssemblyBegin(form, MAT_FINAL_ASSEMBLY);
297 PISM_CHK(ierr, "MatAssemblyBegin");
298
299 ierr = MatAssemblyEnd(form, MAT_FINAL_ASSEMBLY);
300 PISM_CHK(ierr, "MatAssemblyEnd");
301}
302
303} // end of namespace inverse
304} // end of namespace pism
Makes sure that we call begin_access() and end_access() for all accessed array::Arrays.
Definition Array.hh:66
void set(double c)
Result: v[j] <- c for all j.
Definition Array.cc:659
bool grounded_ice(int i, int j) const
Definition CellType.hh:46
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
fem::ElementIterator m_element_index
std::shared_ptr< const Grid > m_grid
virtual void valueAt(array::Scalar &x, double *OUTPUT)
virtual void gradientAt(array::Scalar &x, array::Scalar &gradient)
virtual void dot(array::Scalar &a, array::Scalar &b, double *OUTPUT)
Computes the inner product .
#define PISM_CHK(errcode, name)
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 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