Reference documentation for deal.II version 9.0.0
|
#include <deal.II/base/quadrature_lib.h>
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 ()=default |
Quadrature & | operator= (const Quadrature< dim > &) |
Quadrature & | operator= (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 |
std::conditional< dim==1, std::array< Quadrature< 1 >, dim >, const std::array< Quadrature< 1 >, dim > & >::type | get_tensor_basis () const |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | subscribe (const char *identifier=nullptr) const |
void | unsubscribe (const char *identifier=nullptr) const |
unsigned int | n_subscriptions () 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 > | |
typedef Quadrature< dim-1 > | SubQuadrature |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (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 |
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 727 of file quadrature_lib.h.
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.
radial_quadrature | Base quadrature to use in the radial direction |
angular_quadrature | Base quadrature to use in the angular direction |
Definition at line 1278 of file quadrature_lib.cc.
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.
n | Order of QGauss quadrature |
Definition at line 1306 of file quadrature_lib.cc.