Reference documentation for deal.II version 9.3.3
\(\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 | Protected Attributes | List of all members
QTelles< dim > Class Template Reference

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

Inheritance diagram for QTelles< dim >:
[legend]

Public Types

using SubQuadrature = Quadrature< dim - 1 >
 

Public Member Functions

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

Protected Attributes

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
 

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
 
using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 
std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 
static std::mutex mutex
 
static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 
void check_no_subscribers () const noexcept
 

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.

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

◆ QTelles() [1/3]

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/3]

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.

◆ QTelles() [3/3]

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

Definition at line 883 of file quadrature_lib.cc.

Member Function Documentation

◆ 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
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 using the BOOST serialization library.

◆ 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 >, conststd::array< Quadrature< 1 >, dim > & >::type Quadrature< dim >::get_tensor_basis
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.

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 283 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 289 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 298 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 304 of file quadrature.h.


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