Reference documentation for deal.II version GIT relicensing-446-ge11b936273 2024-04-20 12:30: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
QWitherdenVincentSimplex< dim > Class Template Reference

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

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

Public Member Functions

 QWitherdenVincentSimplex (const unsigned int n_points_1D, const bool use_odd_order=true)
 
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 QWitherdenVincentSimplex< dim >

Witherden-Vincent rules for simplex entities.

Like QGauss, users should specify a number n_points_1d as an indication of what polynomial degree to be integrated exactly (e.g., for \(n\) points, the rule can integrate polynomials of degree \(2 n - 1\) exactly). Additionally, since these rules were derived for simplices, there are also even-ordered rules (i.e., they integrate polynomials of degree \(2 n\)) available which do not have analogous 1d rules.

The given value for n_points_1d = 1, 2, 3, 4, 5, 6, 7 (where the last two are only implemented in 2d) results in the following number of quadrature points in 2d and 3d:

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

These rules match the ones listed for Witherden-Vincent in the quadpy [178] library and were first described in [196].

Note
Some rules (2d 2 odd and 3d 2 even) do not yet exist and instead a higher-order rule is used in their place.

Also see Simplex support.

Definition at line 946 of file quadrature_lib.h.

Constructor & Destructor Documentation

◆ QWitherdenVincentSimplex()

template<int dim>
QWitherdenVincentSimplex< dim >::QWitherdenVincentSimplex ( const unsigned int  n_points_1D,
const bool  use_odd_order = true 
)
explicit

Constructor taking the equivalent number of quadrature points in 1d n_points_1d and boolean indicating whether the rule should be order \(2 n - 1\) or \(2 n\): see the general documentation of this class for more information.

Definition at line 1806 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: