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

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

Inheritance diagram for QTelles< dim >:
[legend]

Public Member Functions

 QTelles (const Quadrature< 1 > &base_quad, const Point< dim > &singularity)
 
 QTelles (const unsigned int n, const Point< dim > &singularity)
 
- 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

template<int dim>
class QTelles< dim >

Telles quadrature of arbitrary order.

The coefficients of these quadrature rules are computed using a non linear change of variables starting from a Gauss-Legendre quadrature formula. This is done using a cubic polynomial, \(n = a x^3 + b x^2 + c x + d\) in order to integrate a singular integral, with singularity at a given point x_0.

We start from a Gauss Quadrature Formula with arbitrary function. Then we apply the cubic variable change. In the paper, J.C.F.Telles:A Self-Adaptive Co-ordinate Transformation For Efficient Numerical Evaluation of General Boundary Element Integrals. International Journal for Numerical Methods in Engineering, vol 24, pages 959–973. year 1987, the author applies the transformation on the reference cell \([-1, 1]\) getting

\begin{align*} n(1) &= 1, \\ n(-1) &= -1, \\ \frac{dn}{dx} &= 0 \text{ at } x = x_0, \\ \frac{d^2n}{dx^2} &= 0 \text{ at } x = x_0 \end{align*}

We get

\begin{align*} a &= \frac{1}{q}, \\ b &= -3 \frac{\bar{\Gamma}}{q}, \\ c &= 3 \frac{\bar{\Gamma}}{q}, \\ d &= -b, \end{align*}

with

\begin{align*} \eta^{*} &= \bar{\eta}^2 - 1, \\ \bar{\Gamma} &= \sqrt[3]{\bar{\eta} \eta^{*} + |\eta^{*} | } + \sqrt[3]{ \bar{\eta} \eta^{*} - |\eta^{*} | } + \bar{\eta}, \\ q &= (\Gamma-\bar{\Gamma})^3 + \bar{\Gamma} \frac{\bar{\Gamma}^2+3}{1+3\bar{\Gamma}^2} \end{align*}

Since the library assumes \([0,1]\) as reference interval, we will map these values on the proper reference interval in the implementation.

This variable change can be used to integrate singular integrals. One example is \(f(x)/|x-x_0|\) on the reference interval \([0,1]\), where \(x_0\) is given at construction time, and is the location of the singularity \(x_0\), and \(f(x)\) is a smooth non singular function.

Singular quadrature formula are rather expensive, nevertheless Telles' quadrature formula are much easier to compute with respect to other singular integration techniques as Lachat-Watson.

We have implemented the case for \(dim = 1\). When we deal the case \(dim >1\) we have computed the quadrature formula has a tensorial product of one dimensional Telles' quadrature formulas considering the different components of the singularity.

The weights and functions for Gauss Legendre formula have been tabulated up to order 12.

Author
Nicola Giuliani, Luca Heltai 2015

Definition at line 467 of file quadrature_lib.h.

Constructor & Destructor Documentation

◆ QTelles() [1/2]

template<int dim>
QTelles< dim >::QTelles ( const Quadrature< 1 > &  base_quad,
const Point< dim > &  singularity 
)

A constructor that takes a quadrature formula and a singular point as argument. The quadrature formula will be mapped using Telles' rule. Make sure that the order of the quadrature rule is appropriate for the singularity in question.

Definition at line 857 of file quadrature_lib.cc.

◆ QTelles() [2/2]

template<int dim>
QTelles< dim >::QTelles ( const unsigned int  n,
const Point< dim > &  singularity 
)

A variant of above constructor that takes as parameters the order n and location of a singularity. A Gauss Legendre quadrature of order n will be used

Definition at line 874 of file quadrature_lib.cc.


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