Reference documentation for deal.II version 9.0.0
|
#include <deal.II/base/quadrature_lib.h>
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) | |
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) |
Static Private Member Functions | |
static unsigned int | quad_size (const Point< dim > singularity, const unsigned int n) |
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 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 295 of file quadrature_lib.h.
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< 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:
|
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.