Reference documentation for deal.II version 9.3.3
\(\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\}}\)
Public Types | Public Member Functions | Protected Attributes | List of all members
QWitherdenVincentSimplex< dim > Class Template Reference

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

Inheritance diagram for QWitherdenVincentSimplex< dim >:
[legend]

Public Types

using SubQuadrature = Quadrature< dim - 1 >
 

Public Member Functions

 QWitherdenVincentSimplex (const unsigned int n_points_1D)
 
Quadrature< dim > compute_affine_transformation (const std::array< Point< dim >, dim+1 > &vertices) const
 
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
 

Protected Attributes

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
 

Subscriptor functionality

Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class.

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
 
using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 
std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 
static std::mutex mutex
 
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 
void check_no_subscribers () const noexcept
 

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). The given value for n_points_1D = 1, 2, 3, 4, 5 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).

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

Definition at line 843 of file quadrature_lib.h.

Member Typedef Documentation

◆ SubQuadrature

template<int dim>
using Quadrature< dim >::SubQuadrature = Quadrature<dim - 1>
inherited

Define an alias for a quadrature that acts on an object of one dimension less. For cells, this would then be a face quadrature.

Definition at line 90 of file quadrature.h.

Constructor & Destructor Documentation

◆ QWitherdenVincentSimplex()

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

Constructor taking the number of quadrature points in 1D direction n_points_1D.

Definition at line 1512 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
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.

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

◆ operator==()

template<int dim>
bool Quadrature< dim >::operator== ( const Quadrature< dim > &  p) const
inherited

Test for equality of two quadratures.

Definition at line 302 of file quadrature.cc.

◆ initialize()

template<int dim>
void Quadrature< dim >::initialize ( const std::vector< Point< dim > > &  points,
const std::vector< double > &  weights 
)
inherited

Set the quadrature points and weights to the values provided in the arguments.

Definition at line 50 of file quadrature.cc.

◆ size()

template<int dim>
unsigned int Quadrature< dim >::size ( ) const
inherited

Number of quadrature points.

◆ point()

template<int dim>
const Point< dim > & Quadrature< dim >::point ( const unsigned int  i) const
inherited

Return the ith quadrature point.

◆ get_points()

template<int dim>
const std::vector< Point< dim > > & Quadrature< dim >::get_points ( ) const
inherited

Return a reference to the whole array of quadrature points.

◆ weight()

template<int dim>
double Quadrature< dim >::weight ( const unsigned int  i) const
inherited

Return the weight of the ith quadrature point.

◆ get_weights()

template<int dim>
const std::vector< double > & Quadrature< dim >::get_weights ( ) const
inherited

Return a reference to the whole array of weights.

◆ memory_consumption()

template<int dim>
std::size_t Quadrature< dim >::memory_consumption
inherited

Determine an estimate for the memory consumption (in bytes) of this object.

Definition at line 311 of file quadrature.cc.

◆ serialize()

template<int dim>
template<class Archive >
void Quadrature< dim >::serialize ( Archive &  ar,
const unsigned int  version 
)
inherited

Write or read the data of this object to or from a stream for the purpose of serialization using the BOOST serialization library.

◆ is_tensor_product()

template<int dim>
bool Quadrature< dim >::is_tensor_product ( ) const
inherited

This function returns true if the quadrature object is a tensor product of one-dimensional formulas and the quadrature points are sorted lexicographically.

◆ get_tensor_basis()

template<int dim>
std::conditional< dim==1, std::array< Quadrature< 1 >, dim >, conststd::array< Quadrature< 1 >, dim > & >::type Quadrature< dim >::get_tensor_basis
inherited

In case the quadrature formula is a tensor product, this function returns the dim one-dimensional basis objects. Otherwise, calling this function is not allowed.

For dim equal to one, we can not return the std::array as a const reference and have to return it by value. In this case, the array will always contain a single element (this).

Note
The actual return type of this function is
std::conditional<dim == 1,
std::array<Quadrature<1>, dim>,
const std::array<Quadrature<1>, dim> &>::type
The type is abbreviated in the online documentation to improve readability of this page.

Definition at line 323 of file quadrature.cc.

Member Data Documentation

◆ quadrature_points

template<int dim>
std::vector<Point<dim> > Quadrature< dim >::quadrature_points
protectedinherited

List of quadrature points. To be filled by the constructors of derived classes.

Definition at line 283 of file quadrature.h.

◆ weights

template<int dim>
std::vector<double> Quadrature< dim >::weights
protectedinherited

List of weights of the quadrature points. To be filled by the constructors of derived classes.

Definition at line 289 of file quadrature.h.

◆ is_tensor_product_flag

template<int dim>
bool Quadrature< dim >::is_tensor_product_flag
protectedinherited

Indicates if this object represents quadrature formula that is a tensor product of one-dimensional formulas. This flag is set if dim==1 or the constructors taking a Quadrature<1> (and possibly a Quadrature<dim-1> object) is called. This implies that the quadrature points are sorted lexicographically.

Definition at line 298 of file quadrature.h.

◆ tensor_basis

template<int dim>
std::unique_ptr<std::array<Quadrature<1>, dim> > Quadrature< dim >::tensor_basis
protectedinherited

Stores the one-dimensional tensor basis objects in case this object can be represented by a tensor product.

Definition at line 304 of file quadrature.h.


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