Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
Public Member Functions | List of all members
QSimplex< dim > Class Template Reference

#include <deal.II/base/quadrature_lib.h>

Inheritance diagram for QSimplex< dim >:
[legend]

Public Member Functions

 QSimplex (const Quadrature< dim > &quad)
 
Quadrature< dim > compute_affine_transformation (const std::array< Point< dim >, dim+1 > &vertices) const
 
- Public Member Functions inherited from Quadrature< dim >
 Quadrature (const unsigned int n_quadrature_points=0)
 
 Quadrature (const SubQuadrature &, const Quadrature< 1 > &)
 
 Quadrature (const Quadrature< dim !=1 ? 1 :0 > &quadrature_1d)
 
 Quadrature (const Quadrature< dim > &q)
 
 Quadrature (Quadrature< dim > &&) noexcept=default
 
 Quadrature (const std::vector< Point< dim >> &points, const std::vector< double > &weights)
 
 Quadrature (const std::vector< Point< dim >> &points)
 
 Quadrature (const Point< dim > &point)
 
virtual ~Quadrature () override=default
 
Quadratureoperator= (const Quadrature< dim > &)
 
Quadratureoperator= (Quadrature< dim > &&)=default
 
bool operator== (const Quadrature< dim > &p) const
 
void initialize (const std::vector< Point< dim >> &points, const std::vector< double > &weights)
 
unsigned int size () const
 
const Point< dim > & point (const unsigned int i) const
 
const std::vector< Point< dim > > & get_points () const
 
double weight (const unsigned int i) const
 
const std::vector< double > & get_weights () const
 
std::size_t memory_consumption () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
bool is_tensor_product () const
 
const std::array< Quadrature< 1 >, dim > & get_tensor_basis () const
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&) noexcept
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&) noexcept
 
void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Additional Inherited Members

- Public Types inherited from Quadrature< dim >
using SubQuadrature = Quadrature< dim - 1 >
 
- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 
- Protected Attributes inherited from Quadrature< dim >
std::vector< Point< dim > > quadrature_points
 
std::vector< double > weights
 
bool is_tensor_product_flag
 
std::unique_ptr< std::array< Quadrature< 1 >, dim > > tensor_basis
 

Detailed Description

template<int dim>
class QSimplex< dim >

Given an arbitrary quadrature formula, return one that chops the quadrature points above the hyper-plane defined by \(\sum_i x_i = 1\). In other words, it extracts those quadrature points from the base formula that satisfy \(\sum_i (\mathbf x_q)_i \le 1+10^{-12}\)."

In general the resulting quadrature is not very useful, unless the quadrature you started from has been constructed specifically to integrate over triangles or tetrahedra. This class only ensures that the resulting quadrature formula only has quadrature points in the reference simplex or on its boundary.

No transformation is applied to the weights, and the weights referring to points that live outside the reference simplex are simply discarded.

The main use of this quadrature formula is not to chop tensor product quadratures. Ideally you should pass to this class a quadrature formula constructed directly using points and weights in the reference simplex, capable of integrating on triangles or tetrahedra.

For finite elements based on quadrilaterals and hexahedra, a QSimplex quadrature formula is not very useful on its own. This class is typically used in conjunction with other classes, like QSplit, to patch the reference element using several QSimplex quadrature formulas.

Such quadrature formulas are useful to integrate functions with singularities at certain points, or functions that present jumps along a co-dimension one surface inside the reference element, like in the extended finite element method (XFEM).

Author
Luca Heltai, 2017.

Definition at line 611 of file quadrature_lib.h.

Constructor & Destructor Documentation

◆ QSimplex()

template<int dim>
QSimplex< dim >::QSimplex ( const Quadrature< dim > &  quad)

Construct a quadrature that only contains the points that are in the lower left reference simplex.

Parameters
[in]quadThe input quadrature.

Definition at line 1198 of file quadrature_lib.cc.

Member Function Documentation

◆ compute_affine_transformation()

template<int dim>
Quadrature< dim > QSimplex< dim >::compute_affine_transformation ( const std::array< Point< dim >, dim+1 > &  vertices) const

Return an affine transformation of this quadrature, that can be used to integrate on the simplex identified by vertices.

Both the quadrature point locations and the weights are transformed, so that you can effectively use the resulting quadrature to integrate on the simplex.

The transformation is defined as

\[ x = v_0 + B \hat x \]

where the matrix \(B\) is given by \(B_{ij} = v[j][i]-v[0][i]\).

The weights are scaled with the absolute value of the determinant of \(B\), that is \(J \dealcoloneq |\text{det}(B)|\). If \(J\) is zero, an empty quadrature is returned. This may happen, in two dimensions, if the three vertices are aligned, or in three dimensions if the four vertices are on the same plane.

Parameters
[in]verticesThe vertices of the simplex you wish to integrate on
Returns
A quadrature object that can be used to integrate on the simplex

Definition at line 1220 of file quadrature_lib.cc.


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