PISM, A Parallel Ice Sheet Model  stable v2.1-1-g6902d5502 committed by Ed Bueler on 2023-12-20 08:38:27 -0800
BlatterTestXZ.cc
Go to the documentation of this file.
1 /* Copyright (C) 2020, 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/BlatterTestXZ.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 BlatterTestXZ::BlatterTestXZ(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 softness (enthalpy and pressure values are irrelevant)
43  m_A = m_flow_law->softness(1e5, 0);
44 
45  // surface elevation parameter (1/m)
46  m_alpha = 4e-8;
47 
48  // surface elevation parameter (m)
49  m_s0 = 2000.0;
50 
51  // ice thickness (m)
52  m_H0 = 1000.0;
53 
54  // sliding parameter
55  m_beta = convert(m_sys, 1.0, "kPa year m-1", "Pa s m-1");
56 
57  m_rho = m_config->get_number("constants.ice.density");
58  m_g = m_config->get_number("constants.standard_gravity");
59 }
60 
62  const int *node_type,
63  const double *ice_bottom,
64  const double *sea_level) {
65  (void) face;
66  (void) node_type;
67  (void) ice_bottom;
68  (void) sea_level;
69 
70  return false;
71 }
72 
73 bool BlatterTestXZ::dirichlet_node(const DMDALocalInfo &info,
75  // use Dirichlet BC at x == -Lx and x == Lx
76  return (I.i == 0 or I.i == info.mx - 1);
77 }
78 
79 Vector2d BlatterTestXZ::u_bc(double x, double y, double z) const {
80  (void) y;
81 
82  return blatter_xz_exact(x, z, m_A, m_rho, m_g, m_s0, m_alpha, m_H0, m_beta);
83 }
84 
86  const double *surface,
87  const double *bed,
88  Vector2d *residual) {
89  (void) surface;
90  (void) bed;
91 
92  // compute x and z coordinates of quadrature points
93  double
94  *x = m_work[0],
95  *z = m_work[1];
96  {
97  double
98  *x_nodal = m_work[2],
99  *z_nodal = m_work[3];
100 
101  for (int n = 0; n < fem::q13d::n_chi; ++n) {
102  x_nodal[n] = element.x(n);
103  z_nodal[n] = element.z(n);
104  }
105 
106  element.evaluate(x_nodal, x);
107  element.evaluate(z_nodal, z);
108  }
109 
110  // loop over all quadrature points
111  for (int q = 0; q < element.n_pts(); ++q) {
112  auto W = element.weight(q) / m_scaling;
113 
114  // loop over all test functions
115  for (int t = 0; t < element.n_chi(); ++t) {
116  const auto &psi = element.chi(q, t);
117 
118  residual[t] += W * psi.val * blatter_xz_source(x[q], z[q], m_A, m_rho, m_g,
119  m_s0, m_alpha, m_H0, m_beta);
120  }
121  }
122 }
123 
124 /*!
125  * Basal contribution to the residual.
126  */
128  const fem::Q1Element3Face &face,
129  const double *tauc_nodal,
130  const double *f_nodal,
131  const Vector2d *u_nodal,
132  Vector2d *residual) {
133 
134  // The basal sliding BC contribution:
135  Blatter::residual_basal(element, face, tauc_nodal, f_nodal, u_nodal, residual);
136 
137  // Additional contribution needed to satisfy the manufactured solution:
138  {
139  // compute x and z coordinates of quadrature points
140  double
141  *x = m_work[0],
142  *z = m_work[1];
143  {
144  double
145  *x_nodal = m_work[2],
146  *z_nodal = m_work[3];
147 
148  for (int n = 0; n < fem::q13d::n_chi; ++n) {
149  x_nodal[n] = element.x(n);
150  z_nodal[n] = element.z(n);
151  }
152 
153  face.evaluate(x_nodal, x);
154  face.evaluate(z_nodal, z);
155  }
156 
157  for (int q = 0; q < face.n_pts(); ++q) {
158  auto W = face.weight(q) / m_scaling;
159 
160  // loop over all test functions
161  for (int t = 0; t < element.n_chi(); ++t) {
162  auto psi = face.chi(q, t);
163 
164  residual[t] += - W * psi * blatter_xz_source_bed(x[q], z[q], m_A, m_rho, m_g,
165  m_s0, m_alpha, m_H0, m_beta);
166  }
167  }
168  }
169 }
170 
171 /*!
172  * Top surface contribution to the residual.
173  */
175  const fem::Q1Element3Face &face,
176  Vector2d *residual) {
177 
178  // compute x and z coordinates of quadrature points
179  double
180  *x = m_work[0],
181  *z = m_work[1];
182  {
183  double
184  *x_nodal = m_work[2],
185  *z_nodal = m_work[3];
186 
187  for (int n = 0; n < fem::q13d::n_chi; ++n) {
188  x_nodal[n] = element.x(n);
189  z_nodal[n] = element.z(n);
190  }
191 
192  face.evaluate(x_nodal, x);
193  face.evaluate(z_nodal, z);
194  }
195 
196  for (int q = 0; q < face.n_pts(); ++q) {
197  auto W = face.weight(q) / m_scaling;
198 
199  auto F = blatter_xz_source_surface(x[q], z[q], m_A, m_rho, m_g,
200  m_s0, m_alpha, m_H0, m_beta);
201 
202  // loop over all test functions
203  for (int t = 0; t < element.n_chi(); ++t) {
204  auto psi = face.chi(q, t);
205 
206  residual[t] += - W * psi * F;
207  }
208  }
209 }
210 
211 } // end of namespace stressbalance
212 } // end of namespace pism
const units::System::Ptr m_sys
unit system used by this component
Definition: Component.hh:160
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
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.
BlatterTestXZ(std::shared_ptr< const Grid > grid, int Mz, int coarsening_factor)
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)
Vector2d u_bc(double x, double y, double z) const
bool marine_boundary(int face, const int *node_type, const double *ice_bottom, const double *sea_level)
void residual_source_term(const fem::Q1Element3 &element, const double *surface, const double *bed, Vector2d *residual)
void residual_surface(const fem::Q1Element3 &element, const fem::Q1Element3Face &face, Vector2d *residual)
double m_A
constant ice hardness
bool dirichlet_node(const DMDALocalInfo &info, const fem::Element3::GlobalIndex &I)
virtual 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)
Definition: residual.cc:163
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
double convert(System::Ptr system, double input, const std::string &spec1, const std::string &spec2)
Convert a quantity from unit1 to unit2.
Definition: Units.cc:70
static double F(double SL, double B, double H, double alpha)
Vector2d blatter_xz_source_surface(double x, double z, double A, double rho, double g, double s_0, double alpha, double H, double beta)
Vector2d blatter_xz_source(double x, double z, double A, double rho, double g, double s_0, double alpha, double H, double beta)
Vector2d blatter_xz_exact(double x, double z, double A, double rho, double g, double s_0, double alpha, double H, double beta)
Vector2d blatter_xz_source_bed(double x, double z, double A, double rho, double g, double s_0, double alpha, double H, double beta)
const int I[]
Definition: ssafd_code.cc:24