Reference documentation for deal.II version GIT relicensing-478-g3275795167 2024-04-24 07:10:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Public Member Functions | List of all members
QGaussSimplex< dim > Class Template Reference

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

Inheritance diagram for QGaussSimplex< dim >:
Inheritance graph
[legend]

Public Member Functions

 QGaussSimplex (const unsigned int n_points_1D)
 
template<int spacedim = dim>
Quadrature< spacedim > compute_affine_transformation (const std::array< Point< spacedim >, dim+1 > &vertices) const
 
template<int spacedim = dim>
Quadrature< spacedim > mapped_quadrature (const std::vector< std::array< Point< spacedim >, dim+1 > > &simplices) const
 

Detailed Description

template<int dim>
class QGaussSimplex< dim >

Integration rule for simplex entities.

Users specify a number n_points_1d as an indication of what polynomial degree to be integrated exactly, similarly to the number of points in a QGauss quadrature object, even though the present quadrature formula is not a tensor product. The given value is translated for n_points_1d=1,2,3,4 to following number of quadrature points for 2d and 3d:

For 1d, the quadrature rule degenerates to a QGauss<1>(n_points_1d).

Note
The quadrature rules implemented by this class come from a variety of sources, but all of them have positive quadrature weights.
Several of the schemes implemented by this class are not symmetric with respect to the vertices - i.e., the locations of the mapped quadrature points depends on the numbering of the cell vertices. If you need rules that are independent of the vertex numbering then use QWitherdenVincentSimplex.

Also see Simplex support.

Definition at line 904 of file quadrature_lib.h.

Constructor & Destructor Documentation

◆ QGaussSimplex()

template<int dim>
QGaussSimplex< dim >::QGaussSimplex ( const unsigned int  n_points_1D)
explicit

Constructor taking the number of quadrature points in 1d direction n_points_1d.

Definition at line 1640 of file quadrature_lib.cc.

Member Function Documentation

◆ compute_affine_transformation()

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

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. The present function works also in the codimension one and codimension two case. For instance, when dim=2 and spacedim=3, we can map the quadrature points so that they live on the physical triangle embedded in the three dimensional space. In such a case, the matrix \(B\) is not square anymore.

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 1473 of file quadrature_lib.cc.

◆ mapped_quadrature()

template<int dim>
template<int spacedim>
Quadrature< spacedim > QSimplex< dim >::mapped_quadrature ( const std::vector< std::array< Point< spacedim >, dim+1 > > &  simplices) const
inherited

Given a collection of simplices, this function creates a global quadrature rule on the area covered by the simplices, by mapping the current quadrature on each simplex. A simplex is identified by its vertices, which are stored into an array of Points. Hence, this function can provide quadrature rules on polygons (or polyhedra), as they can be split into simplices.

Parameters
simplicesA std::vector where each entry is an array of dim+1 points, which identifies the vertices of a simplex.
Returns
A Quadrature object on the cell.

Definition at line 1507 of file quadrature_lib.cc.


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