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

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

Inheritance diagram for QDuffy:
[legend]

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)
 
- Public Member Functions inherited from QSimplex< 2 >
 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

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.

Author
Luca Heltai, Nicola Giuliani, 2017.

Definition at line 726 of file quadrature_lib.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.


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