Reference documentation for deal.II version 9.0.0
|
#include <deal.II/base/quadrature_lib.h>
Public Member Functions | |
QGaussLogR (const unsigned int n, const Point< dim > x0=Point< dim >(), const double alpha=1, const bool factor_out_singular_weight=false) | |
QGaussLogR (QGaussLogR< dim > &&) noexcept=default | |
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) |
Protected Attributes | |
const double | fraction |
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 |
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) |
A class for Gauss quadrature with arbitrary logarithmic weighting function. This formula is used to integrate \(\ln(|x-x_0|/\alpha)\;f(x)\) on the interval \([0,1]\), where \(f\) is a smooth function without singularities, and \(x_0\) and \(\alpha\) are given at construction time, and are the location of the singularity \(x_0\) and an arbitrary scaling factor in the singularity.
You have to make sure that the point \(x_0\) is not one of the Gauss quadrature points of order \(N\), otherwise an exception is thrown, since the quadrature weights cannot be computed correctly.
This quadrature formula is rather expensive, since it uses internally two Gauss quadrature formulas of order n to integrate the nonsingular part of the factor, and two GaussLog quadrature formulas to integrate on the separate segments \([0,x_0]\) and \([x_0,1]\). If the singularity is one of the extremes and the factor alpha is 1, then this quadrature is the same as QGaussLog.
The last argument from the constructor allows you to use this quadrature rule in one of two possible ways:
\[ \int_0^1 g(x) dx = \int_0^1 f(x) \ln\left(\frac{|x-x_0|}{\alpha}\right) dx = \sum_{i=0}^N w_i g(q_i) = \sum_{i=0}^N \bar{w}_i f(q_i) \]
Which one of the two sets of weights is provided, can be selected by the factor_out_singular_weight
parameter. If it is false (the default), then the \(\bar{w}_i\) weights are computed, and you should provide only the smooth function \(f(x)\), since the singularity is included inside the quadrature. If the parameter is set to true, then the singularity is factored out of the quadrature formula, and you should provide a function \(g(x)\), which should at least be similar to \(\ln(|x-x_0|/\alpha)\).
Notice that this quadrature rule is worthless if you try to use it for regular functions once you factored out the singularity.
The weights and functions have been tabulated up to order 12.
Definition at line 241 of file quadrature_lib.h.
QGaussLogR< dim >::QGaussLogR | ( | const unsigned int | n, |
const Point< dim > | x0 = Point< dim >() , |
||
const double | alpha = 1 , |
||
const bool | factor_out_singular_weight = false |
||
) |
The constructor takes four arguments: the order of the Gauss formula on each of the segments \([0,x_0]\) and \([x_0,1]\), the actual location of the singularity, the scale factor inside the logarithmic function and a flag that decides whether the singularity is left inside the quadrature formula or it is factored out, to be included in the integrand.
|
defaultnoexcept |
Move constructor. We cannot rely on the move constructor for Quadrature
, since it does not know about the additional member fraction
of this class.
|
protected |
This is the length of interval \((0,origin)\), or 1 if either of the two extremes have been selected.
Definition at line 267 of file quadrature_lib.h.