deal.II version GIT relicensing-2289-g1e5549a87a 2024-12-21 21:30:00+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
QDuffy Class Reference

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

Inheritance diagram for QDuffy:
Inheritance graph
[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)
 
Quadrature< spacedim > compute_affine_transformation (const std::array< Point< spacedim >, dim+1 > &vertices) const
 
Quadrature< spacedim > mapped_quadrature (const std::vector< std::array< Point< spacedim >, dim+1 > > &simplices) const
 

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 802 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 1570 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 1598 of file quadrature_lib.cc.

Member Function Documentation

◆ compute_affine_transformation()

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

◆ mapped_quadrature()

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


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