Reference documentation for deal.II version Git 11de1224af 2020-11-24 16:17:24 -0500
\(\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 | Static Public Member Functions | Protected Attributes | List of all members
QDuffy Class Reference

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

Inheritance diagram for QDuffy:
[legend]

Public Types

using SubQuadrature = Quadrature< dim - 1 >
 

Public Member Functions

 QDuffy (const Quadrature< 1 > &radial_quadrature, const Quadrature< 1 > &angular_quadrature, const double beta=1.0)
 
 QDuffy (const unsigned int n, const double beta)
 
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
 
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
 

Static Public Member Functions

static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Protected Attributes

std::vector< Point< dim > > quadrature_points
 
std::vector< doubleweights
 
bool is_tensor_product_flag
 
std::unique_ptr< std::array< Quadrature< 1 >, dim > > tensor_basis
 

Detailed Description

A quadrature that implements the Duffy transformation from a square to a triangle to integrate singularities in the origin of the reference simplex.

The Duffy transformation is defined as

\[ \begin{pmatrix} x\\ y \end{pmatrix} = \begin{pmatrix} \hat x^\beta (1-\hat y)\\ \hat x^\beta \hat y \end{pmatrix} \]

with determinant of the Jacobian equal to \(J= \beta \hat x^{2\beta-1}\). Such transformation maps the reference square \([0,1]\times[0,1]\) to the reference simplex, by collapsing the left side of the square and squeezing quadrature points towards the origin, and then shearing the resulting triangle to the reference one. This transformation shows good convergence properties when \(\beta = 1\) with singularities of order \(1/R\) in the origin, but different \(\beta\) values can be selected to increase convergence and/or accuracy when higher order Gauss rules are used (see "Generalized Duffy transformation for integrating vertex singularities", S. E. Mousavi, N. Sukumar, Computational Mechanics 2009).

When \(\beta = 1\), this transformation is also known as the Lachat-Watson transformation.

Definition at line 719 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

◆ QDuffy() [1/2]

QDuffy::QDuffy ( const Quadrature< 1 > &  radial_quadrature,
const Quadrature< 1 > &  angular_quadrature,
const double  beta = 1.0 
)

Constructor that allows the specification of different quadrature rules along the "radial" and "angular" directions.

Since this quadrature is not based on a Polar change of coordinates, it is not fully proper to talk about radial and angular directions. However, the effect of the Duffy transformation is similar to a polar change of coordinates, since the resulting quadrature points are aligned radially with respect to the singularity.

Parameters
radial_quadratureBase quadrature to use in the radial direction
angular_quadratureBase quadrature to use in the angular direction
betaExponent used in the transformation

Definition at line 1284 of file quadrature_lib.cc.

◆ QDuffy() [2/2]

QDuffy::QDuffy ( const unsigned int  n,
const double  beta 
)

Call the above constructor with QGauss<1>(n) quadrature formulas for both the radial and angular quadratures.

Parameters
nOrder of QGauss quadrature
betaExponent used in the transformation

Definition at line 1312 of file quadrature_lib.cc.

Member Function Documentation

◆ compute_affine_transformation()

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 1220 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 ( ) const
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.

◆ 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 >, const std::array< Quadrature< 1 >, dim > & >::type Quadrature< dim >::get_tensor_basis ( ) const
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.

◆ subscribe()

void Subscriptor::subscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Subscribes a user of the object by storing the pointer validity. The subscriber may be identified by text supplied as identifier.

Definition at line 136 of file subscriptor.cc.

◆ unsubscribe()

void Subscriptor::unsubscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Unsubscribes a user from the object.

Note
The identifier and the validity pointer must be the same as the one supplied to subscribe().

Definition at line 156 of file subscriptor.cc.

◆ n_subscriptions()

unsigned int Subscriptor::n_subscriptions ( ) const
inlineinherited

Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.

Definition at line 300 of file subscriptor.h.

◆ list_subscribers() [1/2]

template<typename StreamType >
void Subscriptor::list_subscribers ( StreamType &  stream) const
inlineinherited

List the subscribers to the input stream.

Definition at line 317 of file subscriptor.h.

◆ list_subscribers() [2/2]

void Subscriptor::list_subscribers ( ) const
inherited

List the subscribers to deallog.

Definition at line 204 of file subscriptor.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 282 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 288 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 297 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 303 of file quadrature.h.


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