20 #include "pism/stressbalance/blatter/verification/BlatterTestHalfar.hh"
22 #include "pism/rheology/IsothermalGlen.hh"
24 #include "pism/stressbalance/blatter/verification/manufactured_solutions.hh"
27 namespace stressbalance {
30 int Mz,
int coarsening_factor)
31 :
Blatter(grid, Mz, coarsening_factor) {
53 m_g =
m_config->get_number(
"constants.standard_gravity");
59 return I.i == 0 or
I.i == info.mx - 1;
69 if (fabs(x) >=
m_R0) {
72 return m_H0 * pow(1.0 - pow(fabs(x) /
m_R0, 4.0 / 3.0), 3.0 / 7.0);
76 return u_bc(x, 0.0, z).
u;
80 const double *surface,
97 x_nodal[
n] = element.
x(
n);
98 z_nodal[
n] = element.
z(
n);
106 for (
int q = 0; q < element.
n_pts(); ++q) {
112 for (
int t = 0; t < element.
n_chi(); ++t) {
113 const auto &psi = element.
chi(q, t);
115 residual[t] += W * psi.val *
F;
122 const double *surface_nodal,
123 const double *z_nodal,
124 const double *sl_nodal,
126 (void) surface_nodal;
137 for (
int n = 0;
n < element.
n_chi(); ++
n) {
138 x_nodal[
n] = element.
x(
n);
146 for (
int q = 0; q < face.
n_pts(); ++q) {
158 for (
int t = 0; t < element.
n_chi(); ++t) {
159 auto psi = face.
chi(q, t);
161 residual[t] += W * psi *
F * N;
177 double *x_nodal =
m_work[1];
180 x_nodal[
n] = element.
x(
n);
186 for (
int q = 0; q < face.
n_pts(); ++q) {
192 for (
int t = 0; t < element.
n_chi(); ++t) {
193 auto psi = face.
chi(q, t);
195 residual[t] += - W * psi *
F;
205 const double *tauc_nodal,
206 const double *f_nodal,
216 double *x_nodal =
m_work[1];
219 x_nodal[
n] = element.
x(
n);
225 for (
int q = 0; q < face.
n_pts(); ++q) {
231 for (
int t = 0; t < element.
n_chi(); ++t) {
232 auto psi = face.
chi(q, t);
234 residual[t] += - W * psi *
F;
240 const double *tauc_nodal,
241 const double *f_nodal,
const Config::ConstPtr m_config
configuration database used by this component
const std::shared_ptr< const Grid > m_grid
grid used by this component
This class represents a 2D vector field (such as ice velocity) at a certain grid point.
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.
const Germ & chi(unsigned int q, unsigned int k) const
int n_pts() const
Number of quadrature points.
double weight(unsigned int q) const
Weight of the quadrature point q
double weight(unsigned int q) const
const double & chi(unsigned int q, unsigned int k) const
void evaluate(const T *x, T *result) const
const Vector3 & normal(unsigned int q) const
int n_pts() const
Number of quadrature points.
Isothermal Glen ice allowing extra customization.
double H_exact(double x) const
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]
std::shared_ptr< rheology::FlowLaw > m_flow_law
EnthalpyConverter::Ptr m_EC
const int n_chi
Number of shape functions on a Q1 element.
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)