LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
|
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 | |
QuadRule & | operator= (const QuadRule &)=default |
QuadRule & | operator= (QuadRule &&)=default |
~QuadRule ()=default | |
Private Attributes | |
base::RefEl | ref_el_ |
quadDegree_t | degree_ {0} |
Eigen::MatrixXd | points_ |
Eigen::VectorXd | weights_ |
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.
To create a suitable QuadRule object, use dedicated utility functions:
Definition at line 58 of file quad_rule.h.
|
default |
|
default |
|
default |
|
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.
|
inlineexplicit |
Construct a new quadrature rule by specifying reference element, points, weights and degree of exactness explicitly.
ref_el | The reference element for which the quadrature rule is. |
points | The points of the quadrature rule, a matrix of size ref_el.Dimension() x num_points that contains the points as column vectors |
weights | The weights of the quadrature rule, a vector of length num_points |
degree | The 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_.
|
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:
Definition at line 129 of file quad_rule.h.
References degree_.
Referenced by Order().
|
inline |
Return the total number \( n \) of quadrature points (num of columns of points/weights)
Definition at line 158 of file quad_rule.h.
References points_.
Referenced by lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::Eval(), lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::Eval(), lf::fe::MassElementMatrixProvider< SCALAR, REACTION_COEFF >::Eval(), lf::fe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::Eval(), lf::fe::MassEdgeMatrixProvider< SCALAR, COEFF, EDGESELECTOR >::Eval(), lf::fe::ScalarLoadEdgeVectorProvider< SCALAR, FUNCTOR, EDGESELECTOR >::Eval(), lf::fe::FeHierarchicSegment< SCALAR >::EvaluationNodes(), lf::fe::FeHierarchicTria< SCALAR >::EvaluationNodes(), lf::fe::FeHierarchicSegment< SCALAR >::NodalValuesToDofs(), lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToEdgeDofs(), lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToFaceDofs(), lf::fe::FeHierarchicSegment< SCALAR >::NumEvaluationNodes(), and lf::fe::FeHierarchicTria< SCALAR >::NumEvaluationNodes().
|
inline |
Return the order of the quadrature rule.
The order is the degree of (polynomial) exactness + 1
Definition at line 137 of file quad_rule.h.
References Degree().
|
inline |
All quadrature points \( \begin{pmatrix} \vec{\xi}_0, \ldots, \vec{\xi}_{n-1} \end{pmatrix} \) as column vectors.
RefEl().Dimension() x
\( n \) that contains the point coordinates as column vectors. Definition at line 145 of file quad_rule.h.
References points_.
Referenced by lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::Eval(), lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::Eval(), lf::fe::MassElementMatrixProvider< SCALAR, REACTION_COEFF >::Eval(), lf::fe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::Eval(), lf::fe::MassEdgeMatrixProvider< SCALAR, COEFF, EDGESELECTOR >::Eval(), lf::fe::ScalarLoadEdgeVectorProvider< SCALAR, FUNCTOR, EDGESELECTOR >::Eval(), lf::fe::FeHierarchicSegment< SCALAR >::EvaluationNodes(), lf::fe::FeHierarchicTria< SCALAR >::EvaluationNodes(), lf::fe::FeHierarchicSegment< SCALAR >::NodalValuesToDofs(), and lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToFaceDofs().
|
inline |
Output function controlled by variable out_ctrl;.
o | The stream to write to. |
out_ctrl | Determines 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<<().
|
inline |
The reference element \( \hat{K} \) over which this QuadRule integrates.
Definition at line 104 of file quad_rule.h.
References ref_el_.
Referenced by lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::ReactionDiffusionElementMatrixProvider(), and lf::uscalfe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::ScalarLoadElementVectorProvider().
|
inline |
All quadrature weights \( \begin{pmatrix} \omega_0, \ldots, \omega_{n-1} \end{pmatrix} \) as a vector.
Definition at line 152 of file quad_rule.h.
References weights_.
Referenced by lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::Eval(), lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::Eval(), lf::fe::MassElementMatrixProvider< SCALAR, REACTION_COEFF >::Eval(), lf::fe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::Eval(), lf::fe::MassEdgeMatrixProvider< SCALAR, COEFF, EDGESELECTOR >::Eval(), lf::fe::ScalarLoadEdgeVectorProvider< SCALAR, FUNCTOR, EDGESELECTOR >::Eval(), lf::fe::FeHierarchicSegment< SCALAR >::NodalValuesToDofs(), and lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToFaceDofs().
|
private |
Definition at line 178 of file quad_rule.h.
Referenced by Degree().
|
private |
Definition at line 179 of file quad_rule.h.
Referenced by NumPoints(), Points(), PrintInfo(), and QuadRule().
|
private |
Definition at line 177 of file quad_rule.h.
Referenced by QuadRule(), and RefEl().
|
private |
Definition at line 180 of file quad_rule.h.
Referenced by PrintInfo(), QuadRule(), and Weights().