PISM, A Parallel Ice Sheet Model  stable v2.1-1-g6902d5502 committed by Ed Bueler on 2023-12-20 08:38:27 -0800
BlatterTestvanderVeen.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/BlatterTestvanderVeen.hh"
21 #include "pism/util/node_types.hh"
22 
23 #include "pism/rheology/IsothermalGlen.hh"
24 
25 #include "pism/stressbalance/blatter/verification/manufactured_solutions.hh"
26 
27 namespace pism {
28 namespace stressbalance {
29 
30 BlatterTestvanderVeen::BlatterTestvanderVeen(std::shared_ptr<const Grid> grid,
31  int Mz, int coarsening_factor)
32  : Blatter(grid, Mz, coarsening_factor) {
33 
34  // use the isothermal Glen flow law
35  m_flow_law.reset(new rheology::IsothermalGlen("stress_balance.blatter.", *m_config, m_EC));
36 
37  // make sure we use the same Glen exponent
38  assert(m_flow_law->exponent() == 3.0);
39 
40  // make sure the grid is periodic in the Y direction
41  assert(m_grid->periodicity() == grid::Y_PERIODIC);
42 
43  // store constant ice hardness (enthalpy and pressure values are irrelevant)
44  m_B = m_flow_law->hardness(1e5, 0);
45 
46  // ice thickness (m)
47  m_H0 = 600.0;
48 
49  // ice velocity (m / s)
50  m_V0 = convert(m_sys, 300.0, "m / year", "m / s");
51 
52  m_rho_ice = m_config->get_number("constants.ice.density");
53 
54  // s = alpha * H, where s is surface elevation and H ice thickness
55  double rho_w = m_config->get_number("constants.sea_water.density");
56  m_alpha = 1.0 - m_rho_ice / rho_w;
57 
58  m_g = m_config->get_number("constants.standard_gravity");
59 }
60 
61 bool BlatterTestvanderVeen::dirichlet_node(const DMDALocalInfo &info,
63  (void) info;
64  // use Dirichlet BC at x == 0
65  return I.i == 0;
66 }
67 
69  const int *node_type,
70  const double *ice_bottom,
71  const double *sea_level) {
72  // Ignore sea level and ice bottom elevation:
73  (void) ice_bottom;
74  (void) sea_level;
75 
76  auto nodes = fem::q13d::incident_nodes[face];
77 
78  // number of nodes per face
79  int N = 4;
80 
81  // exclude faces that contain at least one node that is not a part of the boundary
82  for (int n = 0; n < N; ++n) {
83  if (not (node_type[nodes[n]] == NODE_BOUNDARY)) {
84  return false;
85  }
86  }
87 
88  return true;
89 }
90 
91 Vector2d BlatterTestvanderVeen::u_bc(double x, double y, double z) const {
92  (void) y;
93  (void) z;
94 
95  return u_exact(x);
96 }
97 
99  double Q0 = m_H0 * m_V0;
101 }
102 
103 double BlatterTestvanderVeen::H_exact(double x) const {
104  double Q0 = m_H0 * m_V0;
106 }
107 
108 double BlatterTestvanderVeen::b_exact(double x) const {
109  return (m_alpha - 1.0) * H_exact(x);
110 }
111 
112 double BlatterTestvanderVeen::beta_exact(double x) const {
113  double Q0 = m_H0 * m_V0;
115 }
116 
118  const fem::Q1Element3Face &face,
119  const double *surface_nodal,
120  const double *z_nodal,
121  const double *sl_nodal,
122  Vector2d *residual) {
123  (void) surface_nodal;
124  (void) sl_nodal;
125  (void) z_nodal;
126 
127  double Q0 = m_H0 * m_V0;
128 
129  // compute x coordinates of quadrature points
130  double *x = m_work[0];
131  {
132  double *x_nodal = m_work[1];
133 
134  for (int n = 0; n < element.n_chi(); ++n) {
135  x_nodal[n] = element.x(n);
136  }
137 
138  face.evaluate(x_nodal, x);
139  }
140 
141  // loop over all quadrature points
142  for (int q = 0; q < face.n_pts(); ++q) {
143  auto W = face.weight(q) / m_scaling;
144 
145  // the normal direction (1, 0, 0) is hard-wired here
147  m_rho_ice, m_g, m_B);
148 
149  // loop over all test functions
150  for (int t = 0; t < element.n_chi(); ++t) {
151  auto psi = face.chi(q, t);
152 
153  residual[t] += - W * psi * F;
154  }
155  }
156 }
157 
158 /*!
159  * Top surface contribution to the residual.
160  */
162  const fem::Q1Element3Face &face,
163  Vector2d *residual) {
164  double Q0 = m_H0 * m_V0;
165 
166  // compute x coordinates of quadrature points
167  double
168  *x = m_work[0];
169  {
170  double *x_nodal = m_work[1];
171 
172  for (int n = 0; n < fem::q13d::n_chi; ++n) {
173  x_nodal[n] = element.x(n);
174  }
175 
176  face.evaluate(x_nodal, x);
177  }
178 
179  for (int q = 0; q < face.n_pts(); ++q) {
180  auto W = face.weight(q) / m_scaling;
181 
183  m_rho_ice, m_g, m_B);
184 
185  // loop over all test functions
186  for (int t = 0; t < element.n_chi(); ++t) {
187  auto psi = face.chi(q, t);
188 
189  residual[t] += - W * psi * F;
190  }
191  }
192 }
193 
194 } // end of namespace stressbalance
195 } // 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
int n_chi() const
Definition: Element.hh:65
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 x(int n) const
Definition: Element.hh:364
3D Q1 element
Definition: Element.hh:351
Isothermal Glen ice allowing extra customization.
Vector2d u_bc(double x, double y, double z) const
bool dirichlet_node(const DMDALocalInfo &info, const fem::Element3::GlobalIndex &I)
bool marine_boundary(int face, const int *node_type, const double *ice_bottom, const double *sea_level)
BlatterTestvanderVeen(std::shared_ptr< const Grid > grid, int Mz, int coarsening_factor)
void residual_surface(const fem::Q1Element3 &element, const fem::Q1Element3Face &face, 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 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
const unsigned int incident_nodes[n_faces][4]
Definition: FEM.hh:223
@ 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
double blatter_xz_vanderveen_beta(double x, double alpha, double H_0, double Q_0, double rho_i, double g, double B)
double blatter_xz_vanderveen_thickness(double x, double alpha, double H_0, double Q_0, double rho_i, double g, double B)
static double F(double SL, double B, double H, double alpha)
Vector2d blatter_xz_vanderveen_source_lateral(double x, double alpha, double H_0, double Q_0, double rho_i, double g, double B)
@ NODE_BOUNDARY
Definition: node_types.hh:35
Vector2d blatter_xz_vanderveen_source_surface(double x, double alpha, double H_0, double Q_0, double rho_i, double g, double B)
Vector2d blatter_xz_vanderveen_exact(double x, double alpha, double H_0, double Q_0, double rho_i, double g, double B)
const int I[]
Definition: ssafd_code.cc:24