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