Reference documentation for deal.II version 9.0.0
Public Member Functions | List of all members
QSplit< dim > Class Template Reference

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

Inheritance diagram for QSplit< dim >:
[legend]

Public Member Functions

 QSplit (const QSimplex< dim > &base, const Point< dim > &split_point)
 
- 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
 
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
 
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 ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (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)
 

Additional Inherited Members

- Public Types inherited from Quadrature< dim >
typedef Quadrature< dim-1 > SubQuadrature
 
- 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 QSplit< dim >

A quadrature to use when the cell should be split into subregions to integrate using one or more base quadratures.

Author
Luca Heltai, 2017.

Definition at line 764 of file quadrature_lib.h.

Constructor & Destructor Documentation

◆ QSplit()

template<int dim>
QSplit< dim >::QSplit ( const QSimplex< dim > &  base,
const Point< dim > &  split_point 
)

Construct a quadrature formula by splitting the reference hyper cube into the minimum number of simplices that have vertex zero coinciding with split_point, and patch together affine transformations of the base quadrature. The point split_point should be in the reference element, and an exception is thrown if this is not the case.

In two dimensions, the resulting quadrature formula will be composed of two, three, or four triangular quadrature formulas if split_point coincides with one of the vertices, if it lies on one of the edges, or if it is internal to the reference element respectively.

The same is true for the three dimensional case, with six, eight, ten, or twelve tetrahedral quadrature formulas if split_point coincides with one of the vertices, if it lies on one of the edges, on one of the faces, or if it is internal to the reference element respectively.

The resulting quadrature can be used, for example, to integrate functions with integrable singularities at the split point, provided that you select as base quadrature one that can integrate singular points on vertex zero of the reference simplex.

An example usage in dimension two is given by:

const unsigned int order = 5;
QSplit<2> quad(QTrianglePolar(order), Point<2>(.3,.4));

The resulting quadrature will look like the following:

split_quadrature.png
Parameters
baseBase QSimplex quadrature to use
split_pointWhere to split the hyper cube

Definition at line 1313 of file quadrature_lib.cc.


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