PISM, A Parallel Ice Sheet Model  stable v2.1-1-g6902d5502 committed by Ed Bueler on 2023-12-20 08:38:27 -0800
BlatterTestCFBC.hh
Go to the documentation of this file.
1 /* Copyright (C) 2020, 2021, 2022 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 #ifndef BLATTERTESTCFBC_H
21 #define BLATTERTESTCFBC_H
22 
23 #include "pism/stressbalance/blatter/Blatter.hh"
24 
25 namespace pism {
26 namespace stressbalance {
27 
28 /*!
29  * Implements an X-Z flow line setup testing the implementation of lateral (ice margin)
30  * boundary conditions.
31  */
32 class BlatterTestCFBC : public Blatter {
33 public:
34  BlatterTestCFBC(std::shared_ptr<const Grid> grid, int Mz, int coarsening_factor);
35 
36 private:
37  void init_2d_parameters(const Inputs &inputs);
38 
39  bool dirichlet_node(const DMDALocalInfo &info, const fem::Element3::GlobalIndex& I);
40 
41  Vector2d u_bc(double x, double y, double z) const;
42 
43  void residual_source_term(const fem::Q1Element3 &element,
44  const double *surface,
45  const double *bed,
46  Vector2d *residual);
47 
48  void residual_surface(const fem::Q1Element3 &element,
49  const fem::Q1Element3Face &face,
50  Vector2d *residual);
51 
52  void residual_basal(const fem::Q1Element3 &element,
53  const fem::Q1Element3Face &face,
54  const double *tauc_nodal,
55  const double *f_nodal,
56  const Vector2d *u_nodal,
57  Vector2d *residual);
58 
59  void jacobian_basal(const fem::Q1Element3Face &face,
60  const double *tauc_nodal,
61  const double *f_nodal,
62  const Vector2d *u_nodal,
63  double K[2 * fem::q13d::n_chi][2 * fem::q13d::n_chi]);
64  double m_B;
65  double m_g;
66  double m_rho_i;
67  double m_rho_w;
68  double m_L;
69 };
70 
71 } // end of namespace stressbalance
72 } // end of namespace pism
73 
74 #endif /* BLATTERTESTCFBC_H */
std::shared_ptr< const Grid > grid() const
Definition: Component.cc:105
This class represents a 2D vector field (such as ice velocity) at a certain grid point.
Definition: Vector2d.hh:29
3D Q1 element
Definition: Element.hh:351
void residual_source_term(const fem::Q1Element3 &element, const double *surface, const double *bed, Vector2d *residual)
BlatterTestCFBC(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])
void init_2d_parameters(const Inputs &inputs)
Vector2d u_bc(double x, double y, 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)
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)
const int n_chi
Number of shape functions on a Q1 element.
Definition: FEM.hh:213
static double K(double psi_x, double psi_y, double speed, double epsilon)
const int I[]
Definition: ssafd_code.cc:24