PISM, A Parallel Ice Sheet Model 2.2.1-cd005eec8 committed by Constantine Khrulev on 2025-03-07
Loading...
Searching...
No Matches
IP_H1NormFunctional.cc
Go to the documentation of this file.
1// Copyright (C) 2012, 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/IP_H1NormFunctional.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/Scalar.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 array::AccessScope list(x);
41
43
44 // Loop through all LOCAL elements.
45 const int
50
51 for (int j=ys; j<ys+ym; j++) {
52 for (int i=xs; i<xs+xm; i++) {
53 m_element.reset(i, j);
54
55 // Obtain values of x at the quadrature points for the element.
57 if (dirichletBC) {
58 dirichletBC.enforce_homogeneous(m_element, x_e);
59 }
60 m_element.evaluate(x_e, x_q, dxdx_q, dxdy_q);
61
62 for (unsigned int q=0; q<Nq; q++) {
63 auto W = m_element.weight(q);
64 value += W*(m_cL2*x_q[q]*x_q[q]+ m_cH1*(dxdx_q[q]*dxdx_q[q]+dxdy_q[q]*dxdy_q[q]));
65 } // q
66 } // j
67 } // i
68
69 GlobalSum(m_grid->com, &value, OUTPUT, 1);
70}
71
73
74 const unsigned int Nk = fem::q1::n_chi;
75 const unsigned int Nq = m_element.n_pts();
76 const unsigned int Nq_max = fem::MAX_QUADRATURE_SIZE;
77
78 // The value of the objective
79 double value = 0;
80
81 double a_e[Nk];
82 double a_q[Nq_max], dadx_q[Nq_max], dady_q[Nq_max];
83
84 double b_e[Nk];
85 double b_q[Nq_max], dbdx_q[Nq_max], dbdy_q[Nq_max];
86
87 array::AccessScope list{&a, &b};
88
90
91 // Loop through all LOCAL elements.
92 const int
97
98 for (int j=ys; j<ys+ym; j++) {
99 for (int i=xs; i<xs+xm; i++) {
100 m_element.reset(i, j);
101
102 // Obtain values of x at the quadrature points for the element.
103 m_element.nodal_values(a.array(), a_e);
104 if (dirichletBC) {
105 dirichletBC.enforce_homogeneous(m_element, a_e);
106 }
107 m_element.evaluate(a_e, a_q, dadx_q, dady_q);
108
109 m_element.nodal_values(b.array(), b_e);
110 if (dirichletBC) {
111 dirichletBC.enforce_homogeneous(m_element, b_e);
112 }
113 m_element.evaluate(b_e, b_q, dbdx_q, dbdy_q);
114
115 for (unsigned int q=0; q<Nq; q++) {
116 auto W = m_element.weight(q);
117 value += W*(m_cL2*a_q[q]*b_q[q]+ m_cH1*(dadx_q[q]*dbdx_q[q]+dady_q[q]*dbdy_q[q]));
118 } // q
119 } // j
120 } // i
121
122 GlobalSum(m_grid->com, &value, OUTPUT, 1);
123}
124
125
127
128 const unsigned int Nk = fem::q1::n_chi;
129 const unsigned int Nq = m_element.n_pts();
130 const unsigned int Nq_max = fem::MAX_QUADRATURE_SIZE;
131
132 // Clear the gradient before doing anything with it!
133 gradient.set(0);
134
135 double x_e[Nk];
136 double x_q[Nq_max], dxdx_q[Nq_max], dxdy_q[Nq_max];
137
138 double gradient_e[Nk];
139
140 array::AccessScope list{&x, &gradient};
141
143
144 // Loop through all local and ghosted elements.
145 const int
146 xs = m_element_index.xs,
147 xm = m_element_index.xm,
148 ys = m_element_index.ys,
149 ym = m_element_index.ym;
150
151 for (int j=ys; j<ys+ym; j++) {
152 for (int i=xs; i<xs+xm; i++) {
153
154 // Reset the DOF map for this element.
155 m_element.reset(i, j);
156
157 // Obtain values of x at the quadrature points for the element.
158 m_element.nodal_values(x.array(), x_e);
159 if (dirichletBC) {
160 dirichletBC.constrain(m_element);
161 dirichletBC.enforce_homogeneous(m_element, x_e);
162 }
163 m_element.evaluate(x_e, x_q, dxdx_q, dxdy_q);
164
165 // Zero out the element-local residual in prep for updating it.
166 for (unsigned int k=0; k<Nk; k++) {
167 gradient_e[k] = 0;
168 }
169
170 for (unsigned int q=0; q<Nq; q++) {
171 auto W = m_element.weight(q);
172 const double &x_qq=x_q[q];
173 const double &dxdx_qq=dxdx_q[q], &dxdy_qq=dxdy_q[q];
174 for (unsigned int k=0; k<Nk; k++) {
175 gradient_e[k] += 2*W*(m_cL2*x_qq*m_element.chi(q, k).val +
176 m_cH1*(dxdx_qq*m_element.chi(q, k).dx + dxdy_qq*m_element.chi(q, k).dy));
177 } // k
178 } // q
179 m_element.add_contribution(gradient_e, gradient.array());
180 } // j
181 } // i
182}
183
185
186 const unsigned int Nk = fem::q1::n_chi;
187 const unsigned int Nq = m_element.n_pts();
188
189 PetscErrorCode ierr;
190
191 // Zero out the Jacobian in preparation for updating it.
192 ierr = MatZeroEntries(form);
193 PISM_CHK(ierr, "MatZeroEntries");
194
196
197 // Loop through all the elements.
198 const int
199 xs = m_element_index.xs,
200 xm = m_element_index.xm,
201 ys = m_element_index.ys,
202 ym = m_element_index.ym;
203
204 ParallelSection loop(m_grid->com);
205 try {
206 for (int j=ys; j<ys+ym; j++) {
207 for (int i=xs; i<xs+xm; i++) {
208 // Element-local Jacobian matrix (there are Nk vector valued degrees
209 // of freedom per elment, for a total of (2*Nk)*(2*Nk) = 16
210 // entries in the local Jacobian.
211 double K[Nk][Nk];
212
213 // Initialize the map from global to local degrees of freedom for this element.
214 m_element.reset(i, j);
215
216 // Don't update rows/cols where we project to zero.
217 if (zeroLocs) {
218 zeroLocs.constrain(m_element);
219 }
220
221 // Build the element-local Jacobian.
222 ierr = PetscMemzero(K, sizeof(K));
223 PISM_CHK(ierr, "PetscMemzero");
224
225 for (unsigned int q=0; q<Nq; q++) {
226 auto W = m_element.weight(q);
227 for (unsigned int k = 0; k < Nk; k++) { // Test functions
228 const fem::Germ &test_qk = m_element.chi(q, k);
229 for (unsigned int l = 0; l < Nk; l++) { // Trial functions
230 const fem::Germ &test_ql=m_element.chi(q, l);
231 K[k][l] += W*(m_cL2*test_qk.val*test_ql.val +
232 m_cH1*(test_qk.dx*test_ql.dx +
233 test_qk.dy*test_ql.dy));
234 } // l
235 } // k
236 } // q
237 m_element.add_contribution(&K[0][0], form);
238 } // j
239 } // i
240 } catch (...) {
241 loop.failed();
242 }
243 loop.check();
244
245 if (zeroLocs) {
246 zeroLocs.fix_jacobian(form);
247 }
248
249 ierr = MatAssemblyBegin(form, MAT_FINAL_ASSEMBLY);
250 PISM_CHK(ierr, "MatAssemblyBegin");
251
252 ierr = MatAssemblyEnd(form, MAT_FINAL_ASSEMBLY);
253 PISM_CHK(ierr, "MatAssemblyEnd");
254}
255
256} // end of namespace inverse
257} // end of namespace pism
void failed()
Indicates a failure of a parallel section.
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
fem::ElementIterator m_element_index
std::shared_ptr< const Grid > m_grid
virtual void dot(array::Scalar &a, array::Scalar &b, double *OUTPUT)
Computes the inner product .
virtual void gradientAt(array::Scalar &x, array::Scalar &gradient)
virtual void valueAt(array::Scalar &x, double *OUTPUT)
#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