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 | Private Member Functions | Static Private Member Functions | List of all members
QGaussOneOverR< dim > Class Template Reference

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

Inheritance diagram for QGaussOneOverR< dim >:
[legend]

Public Types

using SubQuadrature = Quadrature< dim - 1 >
 

Public Member Functions

 QGaussOneOverR (const unsigned int n, const Point< dim > &singularity, const bool factor_out_singular_weight=false)
 
 QGaussOneOverR (const unsigned int n, const unsigned int vertex_index, const bool factor_out_singular_weight=false)
 
 QGaussOneOverR (const unsigned int n, const Point< 2 > &singularity, const bool factor_out_singularity)
 
 QGaussOneOverR (const unsigned int n, const unsigned int vertex_index, const bool factor_out_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
 

Private Member Functions

unsigned int quad_size (const Point< 2 > &singularity, const unsigned int n)
 

Static Private Member Functions

static unsigned int quad_size (const Point< dim > &singularity, const unsigned int n)
 

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
 
void check_no_subscribers () const noexcept
 
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)
 

Detailed Description

template<int dim>
class QGaussOneOverR< dim >

A class for Gauss quadrature with \(1/R\) weighting function. This formula can be used to integrate \(1/R \ f(x)\) on the reference element \([0,1]^2\), where \(f\) is a smooth function without singularities, and \(R\) is the distance from the point \(x\) to the vertex \(\xi\), given at construction time by specifying its index. Notice that this distance is evaluated in the reference element.

This quadrature formula is obtained from two QGauss quadrature formulas, upon transforming them into polar coordinate system centered at the singularity, and then again into another reference element. This allows for the singularity to be cancelled by part of the Jacobian of the transformation, which contains \(R\). In practice the reference element is transformed into a triangle by collapsing one of the sides adjacent to the singularity. The Jacobian of this transformation contains \(R\), which is removed before scaling the original quadrature, and this process is repeated for the next half element.

Upon construction it is possible to specify whether we want the singularity removed, or not. In other words, this quadrature can be used to integrate \(g(x) = 1/R\ f(x)\), or simply \(f(x)\), with the \(1/R\) factor already included in the quadrature weights.

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

◆ QGaussOneOverR() [1/4]

template<int dim>
QGaussOneOverR< dim >::QGaussOneOverR ( const unsigned int  n,
const Point< dim > &  singularity,
const bool  factor_out_singular_weight = false 
)

This constructor takes three arguments: the order of the Gauss formula, the point of the reference element in which the singularity is located, and whether we include the weighting singular function inside the quadrature, or we leave it in the user function to be integrated.

Traditionally, quadrature formulas include their weighting function, and the last argument is set to false by default. There are cases, however, where this is undesirable (for example when you only know that your singularity has the same order of 1/R, but cannot be written exactly in this way).

In other words, you can use this function in either of the following way, obtaining the same result:

QGaussOneOverR singular_quad(order, q_point, false);
// This will produce the integral of f(x)/R
for(unsigned int i=0; i<singular_quad.size(); ++i)
integral += f(singular_quad.point(i))*singular_quad.weight(i);
// And the same here
QGaussOneOverR singular_quad_noR(order, q_point, true);
// This also will produce the integral of f(x)/R, but 1/R has to
// be specified.
for(unsigned int i=0; i<singular_quad.size(); ++i) {
double R = (singular_quad_noR.point(i)-cell->vertex(vertex_id)).norm();
integral += f(singular_quad_noR.point(i))*singular_quad_noR.weight(i)/R;
}
double norm(const FEValuesBase< dim > &fe, const ArrayView< const std::vector< Tensor< 1, dim > > > &Du)
Definition: divergence.h:472

◆ QGaussOneOverR() [2/4]

template<int dim>
QGaussOneOverR< dim >::QGaussOneOverR ( const unsigned int  n,
const unsigned int  vertex_index,
const bool  factor_out_singular_weight = false 
)

The constructor takes three arguments: the order of the Gauss formula, the index of the vertex where the singularity is located, and whether we include the weighting singular function inside the quadrature, or we leave it in the user function to be integrated. Notice that this is a specialized version of the previous constructor which works only for the vertices of the quadrilateral.

Traditionally, quadrature formulas include their weighting function, and the last argument is set to false by default. There are cases, however, where this is undesirable (for example when you only know that your singularity has the same order of 1/R, but cannot be written exactly in this way).

In other words, you can use this function in either of the following way, obtaining the same result:

QGaussOneOverR singular_quad(order, vertex_id, false);
// This will produce the integral of f(x)/R
for(unsigned int i=0; i<singular_quad.size(); ++i)
integral += f(singular_quad.point(i))*singular_quad.weight(i);
// And the same here
QGaussOneOverR singular_quad_noR(order, vertex_id, true);
// This also will produce the integral of f(x)/R, but 1/R has to
// be specified.
for(unsigned int i=0; i<singular_quad.size(); ++i) {
double R = (singular_quad_noR.point(i)-cell->vertex(vertex_id)).norm();
integral += f(singular_quad_noR.point(i))*singular_quad_noR.weight(i)/R;
}

◆ QGaussOneOverR() [3/4]

QGaussOneOverR< 2 >::QGaussOneOverR ( const unsigned int  n,
const Point< 2 > &  singularity,
const bool  factor_out_singularity 
)

Definition at line 630 of file quadrature_lib.cc.

◆ QGaussOneOverR() [4/4]

QGaussOneOverR< 2 >::QGaussOneOverR ( const unsigned int  n,
const unsigned int  vertex_index,
const bool  factor_out_singularity 
)

Definition at line 679 of file quadrature_lib.cc.

Member Function Documentation

◆ quad_size() [1/2]

template<int dim>
static unsigned int QGaussOneOverR< dim >::quad_size ( const Point< dim > &  singularity,
const unsigned int  n 
)
staticprivate

Given a quadrature point and a degree n, this function returns the size of the singular quadrature rule, considering whether the point is inside the cell, on an edge of the cell, or on a corner of the cell.

◆ quad_size() [2/2]

unsigned int QGaussOneOverR< 2 >::quad_size ( const Point< 2 > &  singularity,
const unsigned int  n 
)
private

Definition at line 609 of file quadrature_lib.cc.

◆ 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 file: