PISM, A Parallel Ice Sheet Model  stable v2.1-1-g6902d5502 committed by Ed Bueler on 2023-12-20 08:38:27 -0800
BlatterTestHalfar.cc
Go to the documentation of this file.
1 /* Copyright (C) 2020, 2021, 2022, 2023 PISM Authors
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 
20 #include "pism/stressbalance/blatter/verification/BlatterTestHalfar.hh"
21 
22 #include "pism/rheology/IsothermalGlen.hh"
23 
24 #include "pism/stressbalance/blatter/verification/manufactured_solutions.hh"
25 
26 namespace pism {
27 namespace stressbalance {
28 
29 BlatterTestHalfar::BlatterTestHalfar(std::shared_ptr<const Grid> grid,
30  int Mz, int coarsening_factor)
31  : Blatter(grid, Mz, coarsening_factor) {
32 
33  // use the isothermal Glen flow law
34  m_flow_law.reset(new rheology::IsothermalGlen("stress_balance.blatter.", *m_config, m_EC));
35 
36  // make sure we use the same Glen exponent
37  assert(m_flow_law->exponent() == 3.0);
38 
39  // make sure the grid is periodic in the Y direction
40  assert(m_grid->periodicity() == grid::Y_PERIODIC);
41 
42  // store constant ice hardness (enthalpy and pressure values are irrelevant)
43  m_B = m_flow_law->hardness(1e5, 0);
44 
45  // dome radius
46  m_R0 = 750e3;
47 
48  // ice thickness (m)
49  m_H0 = 3600.0;
50 
51  m_rho = m_config->get_number("constants.ice.density");
52 
53  m_g = m_config->get_number("constants.standard_gravity");
54 }
55 
56 bool BlatterTestHalfar::dirichlet_node(const DMDALocalInfo &info,
58  // use Dirichlet BC at x == 0 and the "cliff" near the right boundary
59  return I.i == 0 or I.i == info.mx - 1;
60 }
61 
62 Vector2d BlatterTestHalfar::u_bc(double x, double y, double z) const {
63  (void) y;
64 
65  return blatter_xz_halfar_exact(x, z, m_H0, m_R0, m_rho, m_g, m_B);
66 }
67 
68 double BlatterTestHalfar::H_exact(double x) const {
69  if (fabs(x) >= m_R0) {
70  return 0.0;
71  }
72  return m_H0 * pow(1.0 - pow(fabs(x) / m_R0, 4.0 / 3.0), 3.0 / 7.0);
73 }
74 
75 double BlatterTestHalfar::u_exact(double x, double z) const {
76  return u_bc(x, 0.0, z).u;
77 }
78 
80  const double *surface,
81  const double *bed,
82  Vector2d *residual) {
83 
84  (void) surface;
85  (void) bed;
86 
87  // compute x and z coordinates of quadrature points
88  double
89  *x = m_work[0],
90  *z = m_work[1];
91  {
92  double
93  *x_nodal = m_work[2],
94  *z_nodal = m_work[3];
95 
96  for (int n = 0; n < fem::q13d::n_chi; ++n) {
97  x_nodal[n] = element.x(n);
98  z_nodal[n] = element.z(n);
99  }
100 
101  element.evaluate(x_nodal, x);
102  element.evaluate(z_nodal, z);
103  }
104 
105  // loop over all quadrature points
106  for (int q = 0; q < element.n_pts(); ++q) {
107  auto W = element.weight(q) / m_scaling;
108 
109  auto F = blatter_xz_halfar_source(x[q], z[q], m_H0, m_R0, m_rho, m_g, m_B);
110 
111  // loop over all test functions
112  for (int t = 0; t < element.n_chi(); ++t) {
113  const auto &psi = element.chi(q, t);
114 
115  residual[t] += W * psi.val * F;
116  }
117  }
118 }
119 
121  const fem::Q1Element3Face &face,
122  const double *surface_nodal,
123  const double *z_nodal,
124  const double *sl_nodal,
125  Vector2d *residual) {
126  (void) surface_nodal;
127  (void) sl_nodal;
128 
129  // compute x and z coordinates of quadrature points
130  double
131  *x = m_work[0],
132  *z = m_work[1];
133  {
134  double
135  *x_nodal = m_work[2];
136 
137  for (int n = 0; n < element.n_chi(); ++n) {
138  x_nodal[n] = element.x(n);
139  }
140 
141  face.evaluate(x_nodal, x);
142  face.evaluate(z_nodal, z);
143  }
144 
145  // loop over all quadrature points
146  for (int q = 0; q < face.n_pts(); ++q) {
147  auto W = face.weight(q) / m_scaling;
148  auto N3 = face.normal(q);
149  Vector2d N = {N3.x, N3.y};
150 
151  double F = 0.0;
152  if (x[q] > 0.0) {
153  // this lateral BC is not defined at the left (x = 0) boundary
155  }
156 
157  // loop over all test functions
158  for (int t = 0; t < element.n_chi(); ++t) {
159  auto psi = face.chi(q, t);
160 
161  residual[t] += W * psi * F * N;
162  }
163  }
164 }
165 
166 /*!
167  * Top surface contribution to the residual.
168  */
170  const fem::Q1Element3Face &face,
171  Vector2d *residual) {
172 
173  // compute x coordinates of quadrature points
174  double
175  *x = m_work[0];
176  {
177  double *x_nodal = m_work[1];
178 
179  for (int n = 0; n < fem::q13d::n_chi; ++n) {
180  x_nodal[n] = element.x(n);
181  }
182 
183  face.evaluate(x_nodal, x);
184  }
185 
186  for (int q = 0; q < face.n_pts(); ++q) {
187  auto W = face.weight(q) / m_scaling;
188 
190 
191  // loop over all test functions
192  for (int t = 0; t < element.n_chi(); ++t) {
193  auto psi = face.chi(q, t);
194 
195  residual[t] += - W * psi * F;
196  }
197  }
198 }
199 
200 /*!
201  * Basal contribution to the residual.
202  */
204  const fem::Q1Element3Face &face,
205  const double *tauc_nodal,
206  const double *f_nodal,
207  const Vector2d *u_nodal,
208  Vector2d *residual) {
209  (void) tauc_nodal;
210  (void) f_nodal;
211  (void) u_nodal;
212 
213  // compute x coordinates of quadrature points
214  double *x = m_work[0];
215  {
216  double *x_nodal = m_work[1];
217 
218  for (int n = 0; n < fem::q13d::n_chi; ++n) {
219  x_nodal[n] = element.x(n);
220  }
221 
222  face.evaluate(x_nodal, x);
223  }
224 
225  for (int q = 0; q < face.n_pts(); ++q) {
226  auto W = face.weight(q) / m_scaling;
227 
229 
230  // loop over all test functions
231  for (int t = 0; t < element.n_chi(); ++t) {
232  auto psi = face.chi(q, t);
233 
234  residual[t] += - W * psi * F;
235  }
236  }
237 }
238 
240  const double *tauc_nodal,
241  const double *f_nodal,
242  const Vector2d *u_nodal,
243  double K[2 * fem::q13d::n_chi][2 * fem::q13d::n_chi]) {
244  (void) face;
245  (void) tauc_nodal;
246  (void) f_nodal;
247  (void) u_nodal;
248  (void) K;
249  // empty: residual contribution from the basal boundary does not depend on ice velocity
250 }
251 
252 } // end of namespace stressbalance
253 } // end of namespace pism
const Config::ConstPtr m_config
configuration database used by this component
Definition: Component.hh:158
const std::shared_ptr< const Grid > m_grid
grid used by this component
Definition: Component.hh:156
This class represents a 2D vector field (such as ice velocity) at a certain grid point.
Definition: Vector2d.hh:29
void evaluate(const T *x, T *vals, T *dx, T *dy, T *dz) const
Given nodal values, compute the values and partial derivatives at the quadrature points.
Definition: Element.hh:278
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
int n_chi() const
Definition: Element.hh:65
double weight(unsigned int q) const
Weight of the quadrature point q
Definition: Element.hh:85
double weight(unsigned int q) const
Definition: Element.hh:407
const double & chi(unsigned int q, unsigned int k) const
Definition: Element.hh:412
void evaluate(const T *x, T *result) const
Definition: Element.hh:427
const Vector3 & normal(unsigned int q) const
Definition: Element.hh:419
int n_pts() const
Number of quadrature points.
Definition: Element.hh:404
double z(int n) const
Definition: Element.hh:374
double x(int n) const
Definition: Element.hh:364
3D Q1 element
Definition: Element.hh:351
Isothermal Glen ice allowing extra customization.
void residual_basal(const fem::Q1Element3 &element, const fem::Q1Element3Face &face, const double *tauc_nodal, const double *f_nodal, const Vector2d *u_nodal, Vector2d *residual)
void residual_lateral(const fem::Q1Element3 &element, const fem::Q1Element3Face &face, const double *surface_nodal, const double *z_nodal, const double *sl_nodal, Vector2d *residual)
double u_exact(double x, double z) const
bool dirichlet_node(const DMDALocalInfo &info, const fem::Element3::GlobalIndex &I)
void residual_surface(const fem::Q1Element3 &element, const fem::Q1Element3Face &face, Vector2d *residual)
Vector2d u_bc(double x, double y, double z) const
void residual_source_term(const fem::Q1Element3 &element, const double *surface, const double *bed, Vector2d *residual)
BlatterTestHalfar(std::shared_ptr< const Grid > grid, int Mz, int coarsening_factor)
void jacobian_basal(const fem::Q1Element3Face &face, const double *tauc_nodal, const double *f_nodal, const Vector2d *u_nodal, double K[2 *fem::q13d::n_chi][2 *fem::q13d::n_chi])
double m_work[m_n_work][m_Nq]
Definition: Blatter.hh:108
std::shared_ptr< rheology::FlowLaw > m_flow_law
#define n
Definition: exactTestM.c:37
const int n_chi
Number of shape functions on a Q1 element.
Definition: FEM.hh:213
@ Y_PERIODIC
Definition: Grid.hh:51
static double K(double psi_x, double psi_y, double speed, double epsilon)
static double F(double SL, double B, double H, double alpha)
Vector2d blatter_xz_halfar_source_lateral(double x, double z, double H_0, double R_0, double rho_i, double g, double B)
Vector2d blatter_xz_halfar_source(double x, double z, double H_0, double R_0, double rho_i, double g, double B)
Vector2d blatter_xz_halfar_exact(double x, double z, double H_0, double R_0, double rho_i, double g, double B)
Vector2d blatter_xz_halfar_source_surface(double x, double H_0, double R_0, double rho_i, double g, double B)
Vector2d blatter_xz_halfar_source_base(double x, double H_0, double R_0, double rho_i, double g, double B)
const int I[]
Definition: ssafd_code.cc:24