LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
lf::quad::QuadRule Class Reference

Represents a Quadrature Rule over one of the Reference Elements. More...

#include <lf/quad/quad_rule.h>

Public Member Functions

 QuadRule ()
 Default constructor creating an "invalid quadrature rule".
 
 QuadRule (base::RefEl ref_el, Eigen::MatrixXd points, Eigen::VectorXd weights, quadDegree_t degree)
 Construct a new quadrature rule by specifying reference element, points, weights and degree of exactness explicitly.
 
base::RefEl RefEl () const
 The reference element \( \hat{K} \) over which this QuadRule integrates.
 
quadDegree_t Degree () const
 Return the degree of exactness of this Quadrature Rule.
 
unsigned int Order () const
 Return the order of the quadrature rule.
 
const Eigen::MatrixXd & Points () const
 All quadrature points \( \begin{pmatrix} \vec{\xi}_0, \ldots, \vec{\xi}_{n-1} \end{pmatrix} \) as column vectors.
 
const Eigen::VectorXd & Weights () const
 All quadrature weights \( \begin{pmatrix} \omega_0, \ldots, \omega_{n-1} \end{pmatrix} \) as a vector.
 
base::size_type NumPoints () const
 Return the total number \( n \) of quadrature points (num of columns of points/weights)
 
void PrintInfo (std::ostream &o, int out_ctrl=0) const
 Output function controlled by variable out_ctrl;.
 
Default constructors
 QuadRule (const QuadRule &)=default
 
 QuadRule (QuadRule &&)=default
 
QuadRuleoperator= (const QuadRule &)=default
 
QuadRuleoperator= (QuadRule &&)=default
 
 ~QuadRule ()=default
 

Private Attributes

base::RefEl ref_el_
 
quadDegree_t degree_ {0}
 
Eigen::MatrixXd points_
 
Eigen::VectorXd weights_
 

Detailed Description

Represents a Quadrature Rule over one of the Reference Elements.

A Quadrature rule essentially boils down to pairs of quadrature nodes \( \{\vec{\xi}_0, \ldots, \vec{\xi}_{n-1}\} \) and quadrature weights \( \{\omega_1, \ldots, \omega_{n-1}\}\).

The integral of a function \( f \) over the reference element \(K\) is then approximated by

\[ \int_K f(\vec{x}) \, d\vec{x} \approx \sum_{i=0}^{n-1} f(\vec{\xi_i}) \omega_i \]

See Lecture Document Subsection 2.7.5.2 for more explanations.

Use case

template <typename FUNCTOR>
auto localQuadFunction(const lf::mesh::Mesh &mesh,
const lf::quad::QuadRule &quadrule, FUNCTOR &&f) {
// Variable for summing the result
using value_t = std::invoke_result_t<FUNCTOR, Eigen::VectorXd>;
value_t sum_var{};
// Loop over entities of co-dimension 0 = cells
for (const lf::mesh::Entity *entity : mesh.Entities(0)) {
// Check matching of reference element
LF_VERIFY_MSG(entity->RefEl() == quadrule.RefEl(),
"Mismatch of reference element for " << entity->RefEl());
// Obtain geometry information for entity
const lf::geometry::Geometry &geo{*entity->Geometry()};
// Number of quadrature points
const int P = quadrule.NumPoints();
// Quadrature points
const Eigen::MatrixXd zeta_ref{quadrule.Points()};
// Map quadrature points to physical/world coordinates
const Eigen::MatrixXd zeta{geo.Global(zeta_ref)};
// Quadrature weights
const Eigen::VectorXd w_ref{quadrule.Weights()};
// Gramian determinants
const Eigen::VectorXd gram_dets{geo.IntegrationElement(zeta_ref)};
// Iterate over the quadrature points
for (int l = 0; l < P; ++l) {
sum_var += w_ref[l] * f(zeta.col(l)) * gram_dets[l];
}
}
return sum_var;
}

To create a suitable QuadRule object, use dedicated utility functions:

Definition at line 58 of file quad_rule.h.

Constructor & Destructor Documentation

◆ QuadRule() [1/4]

lf::quad::QuadRule::QuadRule ( const QuadRule & )
default

◆ QuadRule() [2/4]

lf::quad::QuadRule::QuadRule ( QuadRule && )
default

◆ ~QuadRule()

lf::quad::QuadRule::~QuadRule ( )
default

◆ QuadRule() [3/4]

lf::quad::QuadRule::QuadRule ( )
inline

Default constructor creating an "invalid quadrature rule".

This default constructor is needed when storing quadrature rules in member variables of classes, whose initialization cannot be done in an initialization section of a constructor.

Definition at line 75 of file quad_rule.h.

◆ QuadRule() [4/4]

lf::quad::QuadRule::QuadRule ( base::RefEl ref_el,
Eigen::MatrixXd points,
Eigen::VectorXd weights,
quadDegree_t degree )
inlineexplicit

Construct a new quadrature rule by specifying reference element, points, weights and degree of exactness explicitly.

Parameters
ref_elThe reference element for which the quadrature rule is.
pointsThe points of the quadrature rule, a matrix of size ref_el.Dimension() x num_points that contains the points as column vectors
weightsThe weights of the quadrature rule, a vector of length num_points
degreeThe degree of exactness of the quadrature rule, see Degree() for more info.

Definition at line 90 of file quad_rule.h.

References lf::base::RefEl::Dimension(), points_, ref_el_, and weights_.

Member Function Documentation

◆ Degree()

quadDegree_t lf::quad::QuadRule::Degree ( ) const
inline

Return the degree of exactness of this Quadrature Rule.

The degree of exactness a quadrature rule is defined as the largest integer \(k\) such that

\[ \forall p \in \mathbb{P}_k, \quad \int_K p(\vec{x}) \, d\vec{x} = \sum_{i=0}^{n-1} \omega_i p(\vec{\xi_i}) \]

here the polynomial space \(\mathbb{P}_k\) is defined differently for every reference element:

  • For base::RefEl::kSegment(), \( \mathbb{P}_k := \mathrm{span} \{ x^a \, | \, 0 \leq a \leq k \} \).
  • For base::RefEl::kTria(), \( \mathbb{P}_k := \mathrm{span} \{ x^a y^b \, | \, 0 \leq a + b \leq k \} \)
  • For base::RefEl::kQuad(), \( \mathbb{P}_k := \mathrm{span} \{ x^a y^b \, | \, 0 \leq \mathrm{max}(a,b) \leq k \} \)
Note
the order of a quadrature rule is the degree of exactness + 1, because it predicts the order of algebraic covergence of the quadrature error when the rule is applied to a smooth integrand.

Definition at line 129 of file quad_rule.h.

References degree_.

Referenced by Order().

◆ NumPoints()

base::size_type lf::quad::QuadRule::NumPoints ( ) const
inline

◆ operator=() [1/2]

QuadRule & lf::quad::QuadRule::operator= ( const QuadRule & )
default

◆ operator=() [2/2]

QuadRule & lf::quad::QuadRule::operator= ( QuadRule && )
default

◆ Order()

unsigned int lf::quad::QuadRule::Order ( ) const
inline

Return the order of the quadrature rule.

The order is the degree of (polynomial) exactness + 1

See also
Degree()

Definition at line 137 of file quad_rule.h.

References Degree().

◆ Points()

const Eigen::MatrixXd & lf::quad::QuadRule::Points ( ) const
inline

◆ PrintInfo()

void lf::quad::QuadRule::PrintInfo ( std::ostream & o,
int out_ctrl = 0 ) const
inline

Output function controlled by variable out_ctrl;.

Parameters
oThe stream to write to.
out_ctrlDetermines the level of detail of the printed info (see below).

If out_ctrl > 0 weights and nodes are printed, otherwise just the number of nodes are output.

Definition at line 169 of file quad_rule.h.

References points_, and weights_.

Referenced by lf::quad::operator<<().

◆ RefEl()

base::RefEl lf::quad::QuadRule::RefEl ( ) const
inline

◆ Weights()

const Eigen::VectorXd & lf::quad::QuadRule::Weights ( ) const
inline

Member Data Documentation

◆ degree_

quadDegree_t lf::quad::QuadRule::degree_ {0}
private

Definition at line 178 of file quad_rule.h.

Referenced by Degree().

◆ points_

Eigen::MatrixXd lf::quad::QuadRule::points_
private

Definition at line 179 of file quad_rule.h.

Referenced by NumPoints(), Points(), PrintInfo(), and QuadRule().

◆ ref_el_

base::RefEl lf::quad::QuadRule::ref_el_
private

Definition at line 177 of file quad_rule.h.

Referenced by QuadRule(), and RefEl().

◆ weights_

Eigen::VectorXd lf::quad::QuadRule::weights_
private

Definition at line 180 of file quad_rule.h.

Referenced by PrintInfo(), QuadRule(), and Weights().


The documentation for this class was generated from the following file: