27#include "pism/util/fem/FEM.hh"
28#include "pism/util/Vector2d.hh"
35namespace array {
class Scalar; }
73 const Germ&
chi(
unsigned int q,
unsigned int k)
const {
73 const Germ&
chi(
unsigned int q,
unsigned int k)
const {
…}
85 double weight(
unsigned int q)
const {
95 for (
unsigned int q = 0; q <
m_Nq; q++) {
105 Element(
const DMDALocalInfo &grid_info,
int Nq,
int n_chi,
int block_size);
121 const std::vector<QuadPoint>& points,
122 const std::vector<double>& W);
168 void reset(
int i,
int j);
194 template <
typename T>
196 for (
unsigned int q = 0; q <
m_Nq; q++) {
202#ifdef __clang_analyzer__
208 vals[q] += psi.
val * x[
k];
209 dx[q] += psi.
dx * x[
k];
210 dy[q] += psi.
dy * x[
k];
227 result[
k] = x_global[j][i];
247 y_global[j][i] += local[
k];
256 Element2(
const DMDALocalInfo &grid_info,
int Nq,
int n_chi,
int block_size);
268 double dx,
double dy,
284 template <
typename T>
285 void evaluate(
const T *x, T *vals, T *dx, T *dy, T *dz)
const {
286 for (
unsigned int q = 0; q <
m_Nq; q++) {
293 vals[q] += psi.
val * x[
k];
294 dx[q] += psi.
dx * x[
k];
295 dy[q] += psi.
dy * x[
k];
296 dz[q] += psi.
dz * x[
k];
285 void evaluate(
const T *x, T *vals, T *dx, T *dy, T *dz)
const {
…}
324 result[
n] = x_global[I.j][I.i][I.k];
342 y_global[I.j][I.i][I.k] += local[
n];
346 Element3(
const DMDALocalInfo &grid_info,
int Nq,
int n_chi,
int block_size);
368 void reset(
int i,
int j,
int k,
const double *
z);
371 double x(
int n)
const {
376 double y(
int n)
const {
381 double z(
int n)
const {
419 const double&
chi(
unsigned int q,
unsigned int k)
const {
419 const double&
chi(
unsigned int q,
unsigned int k)
const {
…}
431 void reset(
int face,
const double *z);
433 template <
typename T>
435 for (
unsigned int q = 0; q <
m_Nq; q++) {
Describes the PISM grid and the distribution of data across processors.
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)
Given nodal values, compute the values and partial derivatives at the quadrature points.
void local_to_global(unsigned int k, int &i, int &j) const
Convert a local degree of freedom index k to a global degree of freedom index (i,j).
void reset(int i, int j)
Initialize the Element to element (i, j) for the purposes of inserting into global residual and Jacob...
std::vector< double > m_side_lengths
double side_length(unsigned int side) const
std::vector< Vector2d > m_normals
void add_contribution(const T *local, T **y_global) const
Add the values of element-local contributions y to the global vector y_global.
void nodal_values(T const *const *x_global, T *result) const
Extract nodal values for the element (i,j) from global array x_global into the element-local array re...
void nodal_values(const array::Scalar &x_global, int *result) const
Get nodal values of an integer mask.
unsigned int n_sides() const
Vector2d normal(unsigned int side) const
GlobalIndex local_to_global(unsigned int n) const
void nodal_values(T const *const *const *x_global, T *result) const
Extract nodal values for the element (i,j,k) from global array x_global into the element-local array ...
std::vector< int > m_k_offset
void add_contribution(const T *local, T ***y_global) const
Add the values of element-local contributions y to the global vector y_global.
GlobalIndex local_to_global(int i, int j, int k, unsigned int n) const
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.
void evaluate(const T *x, T *result) const
Given nodal values, compute the values at quadrature points.
std::vector< int > m_i_offset
void initialize(const double J[3][3], ShapeFunction f, unsigned int n_chi, const std::vector< QuadPoint > &points, const std::vector< double > &W)
Initialize shape function values and quadrature weights of a 2D physical element.
std::vector< double > m_weights
Quadrature weights for a particular physical element.
Germ(* ShapeFunction)(unsigned int k, const QuadPoint &p)
double weight(unsigned int q) const
Weight of the quadrature point q
std::vector< MatStencil > m_col
const Germ & chi(unsigned int q, unsigned int k) const
const unsigned int m_n_chi
void add_contribution(const double *K, Mat J) const
Add Jacobian contributions.
std::vector< Germ > m_germs
void mark_row_invalid(unsigned int k)
Mark that the row corresponding to local degree of freedom k should not be updated when inserting int...
const unsigned int m_Nq
Number of quadrature points.
unsigned int n_pts() const
Number of quadrature points.
static const int m_invalid_dof
std::vector< MatStencil > m_row
Stencils for updating entries of the Jacobian matrix.
void mark_col_invalid(unsigned int k)
Mark that the column corresponding to local degree of freedom k should not be updated when inserting ...
std::vector< int > m_j_offset
int m_i
Indices of the current element.
The mapping from global to local degrees of freedom.
P1 element embedded in Q1Element2.
Q1 element with sides parallel to X and Y axes.
std::vector< Vector3 > m_normals
double weight(unsigned int q) const
unsigned int n_pts() const
Number of quadrature points.
const double & chi(unsigned int q, unsigned int k) const
std::vector< double > m_w
std::vector< double > m_chi
void evaluate(const T *x, T *result) const
std::vector< QuadPoint > m_points
const Vector3 & normal(unsigned int q) const
std::vector< double > m_weights
const unsigned int m_n_chi
void reset(int face, const double *z)
void reset(int i, int j, int k, const double *z)
std::vector< QuadPoint > m_points
double m_z_nodal[q13d::n_chi]
std::vector< double > m_w
std::vector< Germ > m_chi
Numerical integration of finite element functions.
const int n_chi
Number of shape functions on a Q1 element.
double dy
Function derivative with respect to y.
double val
Function value.
double dz
Function derivative with respect to z.
double dx
Function deriviative with respect to x.
Struct for gathering the value and derivative of a function at a point.