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
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
42
44
45 // Loop through all LOCAL elements.
46 const int
51
52 for (int j=ys; j<ys+ym; j++) {
53 for (int i=xs; i<xs+xm; i++) {
54 bool all_grounded_ice = (m_ice_mask.grounded_ice(i, j) and
55 m_ice_mask.grounded_ice(i+1, j) and
56 m_ice_mask.grounded_ice(i, j+1) and
57 m_ice_mask.grounded_ice(i+1, j+1));
58
59 if (! all_grounded_ice) {
60 continue;
61 }
62
63 m_element.reset(i, j);
64
65 // Obtain values of x at the quadrature points for the element.
67 if (dirichletBC) {
68 dirichletBC.enforce_homogeneous(m_element, x_e);
69 }
70 m_element.evaluate(x_e, x_q, dxdx_q, dxdy_q);
71
72 for (unsigned int q=0; q<Nq; q++) {
73 auto W = m_element.weight(q);
74 value += W*(m_cL2*x_q[q]*x_q[q]+ m_cH1*(dxdx_q[q]*dxdx_q[q]+dxdy_q[q]*dxdy_q[q]));
75 } // q
76 } // j
77 } // i
78
79 GlobalSum(m_grid->com, &value, OUTPUT, 1);
80}
81
83
84 const unsigned int Nk = fem::q1::n_chi;
85 const unsigned int Nq = m_element.n_pts();
86 const unsigned int Nq_max = fem::MAX_QUADRATURE_SIZE;
87
88 // The value of the objective
89 double value = 0;
90
91 double a_e[Nk];
92 double a_q[Nq_max], dadx_q[Nq_max], dady_q[Nq_max];
93
94 array::AccessScope list{&a, &b, &m_ice_mask};
95
96 double b_e[Nk];
97 double b_q[Nq_max], dbdx_q[Nq_max], dbdy_q[Nq_max];
98
100
101 // Loop through all LOCAL elements.
102 const int
103 xs = m_element_index.lxs,
104 xm = m_element_index.lxm,
105 ys = m_element_index.lys,
106 ym = m_element_index.lym;
107
108 for (int j=ys; j<ys+ym; j++) {
109 for (int i=xs; i<xs+xm; i++) {
110 bool all_grounded_ice = (m_ice_mask.grounded_ice(i, j) and
111 m_ice_mask.grounded_ice(i+1, j) and
112 m_ice_mask.grounded_ice(i, j+1) and
113 m_ice_mask.grounded_ice(i+1, j+1));
114
115 if (! all_grounded_ice) {
116 continue;
117 }
118
119 m_element.reset(i, j);
120
121 // Obtain values of x at the quadrature points for the element.
122 m_element.nodal_values(a.array(), a_e);
123 if (dirichletBC) {
124 dirichletBC.enforce_homogeneous(m_element, a_e);
125 }
126 m_element.evaluate(a_e, a_q, dadx_q, dady_q);
127
128 m_element.nodal_values(b.array(), b_e);
129 if (dirichletBC) {
130 dirichletBC.enforce_homogeneous(m_element, b_e);
131 }
132 m_element.evaluate(b_e, b_q, dbdx_q, dbdy_q);
133
134 for (unsigned int q=0; q<Nq; q++) {
135 auto W = m_element.weight(q);
136 value += W*(m_cL2*a_q[q]*b_q[q]+ m_cH1*(dadx_q[q]*dbdx_q[q]+dady_q[q]*dbdy_q[q]));
137 } // q
138 } // j
139 } // i
140
141 GlobalSum(m_grid->com, &value, OUTPUT, 1);
142}
143
144
146
147 const unsigned int Nk = fem::q1::n_chi;
148 const unsigned int Nq = m_element.n_pts();
149 const unsigned int Nq_max = fem::MAX_QUADRATURE_SIZE;
150
151 // Clear the gradient before doing anything with it!
152 gradient.set(0);
153
154 double x_e[Nk];
155 double x_q[Nq_max], dxdx_q[Nq_max], dxdy_q[Nq_max];
156
157 array::AccessScope list{&x, &gradient, &m_ice_mask};
158
159 double gradient_e[Nk];
160
162
163 // Loop through all local and ghosted elements.
164 const int
165 xs = m_element_index.xs,
166 xm = m_element_index.xm,
167 ys = m_element_index.ys,
168 ym = m_element_index.ym;
169
170 for (int j=ys; j<ys+ym; j++) {
171 for (int i=xs; i<xs+xm; i++) {
172 bool all_grounded_ice = (m_ice_mask.grounded_ice(i, j) and
173 m_ice_mask.grounded_ice(i+1, j) and
174 m_ice_mask.grounded_ice(i, j+1) and
175 m_ice_mask.grounded_ice(i+1, j+1));
176
177 if (! all_grounded_ice) {
178 continue;
179 }
180
181 // Reset the DOF map for this element.
182 m_element.reset(i, j);
183
184 // Obtain values of x at the quadrature points for the element.
185 m_element.nodal_values(x.array(), x_e);
186 if (dirichletBC) {
187 dirichletBC.constrain(m_element);
188 dirichletBC.enforce_homogeneous(m_element, x_e);
189 }
190 m_element.evaluate(x_e, x_q, dxdx_q, dxdy_q);
191
192 // Zero out the element-local residual in prep for updating it.
193 for (unsigned int k=0; k<Nk; k++) {
194 gradient_e[k] = 0;
195 }
196
197 for (unsigned int q=0; q<Nq; q++) {
198 auto W = m_element.weight(q);
199 const double &x_qq=x_q[q];
200 const double &dxdx_qq=dxdx_q[q], &dxdy_qq=dxdy_q[q];
201 for (unsigned int k=0; k<Nk; k++) {
202 gradient_e[k] += 2*W*(m_cL2*x_qq*m_element.chi(q, k).val +
203 m_cH1*(dxdx_qq*m_element.chi(q, k).dx + dxdy_qq*m_element.chi(q, k).dy));
204 } // k
205 } // q
206 m_element.add_contribution(gradient_e, gradient.array());
207 } // j
208 } // i
209}
210
212
213 const unsigned int Nk = fem::q1::n_chi;
214 const unsigned int Nq = m_element.n_pts();
215
216 PetscErrorCode ierr;
217 int i, j;
218
219 // Zero out the Jacobian in preparation for updating it.
220 ierr = MatZeroEntries(form);
221 PISM_CHK(ierr, "MatZeroEntries");
222
224
226
227 // Loop through all the elements.
228 const int
229 xs = m_element_index.xs,
230 xm = m_element_index.xm,
231 ys = m_element_index.ys,
232 ym = m_element_index.ym;
233
234 for (j=ys; j<ys+ym; j++) {
235 for (i=xs; i<xs+xm; i++) {
236 bool all_grounded_ice = (m_ice_mask.grounded_ice(i, j) and
237 m_ice_mask.grounded_ice(i+1, j) and
238 m_ice_mask.grounded_ice(i, j+1) and
239 m_ice_mask.grounded_ice(i+1, j+1));
240 if (! all_grounded_ice) {
241 continue;
242 }
243
244 // Element-local Jacobian matrix (there are Nk vector valued degrees
245 // of freedom per elment, for a total of (2*Nk)*(2*Nk) = 16
246 // entries in the local Jacobian.
247 double K[Nk][Nk];
248
249 // Initialize the map from global to local degrees of freedom for this element.
250 m_element.reset(i, j);
251
252 // Don't update rows/cols where we project to zero.
253 if (zeroLocs) {
254 zeroLocs.constrain(m_element);
255 }
256
257 // Build the element-local Jacobian.
258 ierr = PetscMemzero(K, sizeof(K));
259 PISM_CHK(ierr, "PetscMemzero");
260
261 for (unsigned int q=0; q<Nq; q++) {
262 auto W = m_element.weight(q);
263 for (unsigned int k = 0; k < Nk; k++) { // Test functions
264 const fem::Germ &test_qk=m_element.chi(q, k);
265 for (unsigned int l = 0; l < Nk; l++) { // Trial functions
266 const fem::Germ &test_ql=m_element.chi(q, l);
267 K[k][l] += W*(m_cL2*test_qk.val*test_ql.val
268 + m_cH1*(test_qk.dx*test_ql.dx + test_qk.dy*test_ql.dy));
269 } // l
270 } // k
271 } // q
272 m_element.add_contribution(&K[0][0], form);
273 } // j
274 } // i
275
276 if (zeroLocs) {
277 zeroLocs.fix_jacobian(form);
278 }
279
280 ierr = MatAssemblyBegin(form, MAT_FINAL_ASSEMBLY);
281 PISM_CHK(ierr, "MatAssemblyBegin");
282
283 ierr = MatAssemblyEnd(form, MAT_FINAL_ASSEMBLY);
284 PISM_CHK(ierr, "MatAssemblyEnd");
285}
286
287} // end of namespace inverse
288} // 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
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