Reference documentation for deal.II version 8.5.1
Public Types | Public Member Functions | Protected Attributes | List of all members
Quadrature< dim > Class Template Reference

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

Inheritance diagram for Quadrature< dim >:
[legend]

Public Types

typedef Quadrature< dim-1 > SubQuadrature
 

Public Member Functions

 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 > &&)=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 ()
 
Quadratureoperator= (const Quadrature< dim > &)
 
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)
 
- Public Member Functions inherited from Subscriptor
 Subscriptor ()
 
 Subscriptor (const Subscriptor &)
 
 Subscriptor (Subscriptor &&)
 
virtual ~Subscriptor ()
 
Subscriptoroperator= (const Subscriptor &)
 
Subscriptoroperator= (Subscriptor &&)
 
void subscribe (const char *identifier=0) const
 
void unsubscribe (const char *identifier=0) const
 
unsigned int n_subscriptions () const
 
void list_subscribers () const
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 

Protected Attributes

std::vector< Point< dim > > quadrature_points
 
std::vector< double > weights
 

Additional Inherited Members

- Static Public Member Functions inherited from Subscriptor
static ::ExceptionBaseExcInUse (int arg1, char *arg2, std::string &arg3)
 
static ::ExceptionBaseExcNoSubscriber (char *arg1, char *arg2)
 

Detailed Description

template<int dim>
class Quadrature< dim >

Base class for quadrature formulae in arbitrary dimensions. This class stores quadrature points and weights on the unit line [0,1], unit square [0,1]x[0,1], etc.

There are a number of derived classes, denoting concrete integration formulae. Their names names prefixed by Q. Refer to the list of derived classes for more details.

The schemes for higher dimensions are typically tensor products of the one- dimensional formulae, but refer to the section on implementation detail below.

In order to allow for dimension independent programming, a quadrature formula of dimension zero exists. Since an integral over zero dimensions is the evaluation at a single point, any constructor of such a formula initializes to a single quadrature point with weight one. Access to the weight is possible, while access to the quadrature point is not permitted, since a Point of dimension zero contains no information. The main purpose of these formulae is their use in QProjector, which will create a useful formula of dimension one out of them.

Mathematical background

For each quadrature formula we denote by m, the maximal degree of polynomials integrated exactly. This number is given in the documentation of each formula. The order of the integration error is m+1, that is, the error is the size of the cell to the m+1 by the Bramble- Hilbert Lemma. The number m is to be found in the documentation of each concrete formula. For the optimal formulae QGauss we have \(m = 2N-1\), where N is the constructor parameter to QGauss. The tensor product formulae are exact on tensor product polynomials of degree m in each space direction, but they are still only of m+1st order.

Implementation details

Most integration formulae in more than one space dimension are tensor products of quadrature formulae in one space dimension, or more generally the tensor product of a formula in (dim-1) dimensions and one in one dimension. There is a special constructor to generate a quadrature formula from two others. For example, the QGauss<dim> formulae include Ndim quadrature points in dim dimensions, where N is the constructor parameter of QGauss.

Note
Instantiations for this template are provided for dimensions 0, 1, 2, and 3 (see the section on Template instantiations).
Author
Wolfgang Bangerth, Guido Kanschat, 1998, 1999, 2000, 2005, 2009

Definition at line 81 of file quadrature.h.

Member Typedef Documentation

◆ SubQuadrature

template<int dim>
typedef Quadrature<dim-1> Quadrature< dim >::SubQuadrature

Define a typedef for a quadrature that acts on an object of one dimension less. For cells, this would then be a face quadrature.

Definition at line 88 of file quadrature.h.

Constructor & Destructor Documentation

◆ Quadrature() [1/8]

template<int dim>
Quadrature< dim >::Quadrature ( const unsigned int  n_quadrature_points = 0)
explicit

Constructor.

This constructor is marked as explicit to avoid involuntary accidents like in hp::QCollection<dim> q_collection(3) where hp::QCollection<dim> q_collection(QGauss<dim>(3)) was meant.

Definition at line 45 of file quadrature.cc.

◆ Quadrature() [2/8]

template<int dim>
Quadrature< dim >::Quadrature ( const SubQuadrature< dim > &  q1,
const Quadrature< 1 > &  q2 
)

Build this quadrature formula as the tensor product of a formula in a dimension one less than the present and a formula in one dimension.

SubQuadrature<dim>::type expands to Quadrature<dim-1>.

Definition at line 108 of file quadrature.cc.

◆ Quadrature() [3/8]

template<int dim>
Quadrature< dim >::Quadrature ( const Quadrature< dim !=1 ? 1 :0 > &  quadrature_1d)
explicit

Build this quadrature formula as the dim-fold tensor product of a formula in one dimension.

Assuming that the points in the one-dimensional rule are in ascending order, the points of the resulting rule are ordered lexicographically with x running fastest.

In order to avoid a conflict with the copy constructor in 1d, we let the argument be a 0d quadrature formula for dim==1, and a 1d quadrature formula for all other space dimensions.

Definition at line 208 of file quadrature.cc.

◆ Quadrature() [4/8]

template<int dim>
Quadrature< dim >::Quadrature ( const Quadrature< dim > &  q)

Copy constructor.

Definition at line 242 of file quadrature.cc.

◆ Quadrature() [5/8]

template<int dim>
Quadrature< dim >::Quadrature ( Quadrature< dim > &&  )
default

Move constructor. Construct a new quadrature object by transferring the internal data of another quadrature object.

Note
this constructor is only available if deal.II is configured with C++11 support.

◆ Quadrature() [6/8]

template<int dim>
Quadrature< dim >::Quadrature ( const std::vector< Point< dim > > &  points,
const std::vector< double > &  weights 
)

Construct a quadrature formula from given vectors of quadrature points (which should really be in the unit cell) and the corresponding weights. You will want to have the weights sum up to one, but this is not checked.

Definition at line 66 of file quadrature.cc.

◆ Quadrature() [7/8]

template<int dim>
Quadrature< dim >::Quadrature ( const std::vector< Point< dim > > &  points)

Construct a dummy quadrature formula from a list of points, with weights set to infinity. The resulting object is therefore not meant to actually perform integrations, but rather to be used with FEValues objects in order to find the position of some points (the quadrature points in this object) on the transformed cell in real space.

Definition at line 79 of file quadrature.cc.

◆ Quadrature() [8/8]

template<int dim>
Quadrature< dim >::Quadrature ( const Point< dim > &  point)

Constructor for a one-point quadrature. Sets the weight of this point to one.

Definition at line 91 of file quadrature.cc.

◆ ~Quadrature()

template<int dim>
Quadrature< dim >::~Quadrature ( )
virtual

Virtual destructor.

Definition at line 273 of file quadrature.cc.

Member Function Documentation

◆ operator=()

template<int dim>
Quadrature< dim > & Quadrature< dim >::operator= ( const Quadrature< dim > &  q)

Assignment operator. Copies contents of weights and quadrature_points as well as size.

Definition at line 252 of file quadrature.cc.

◆ operator==()

template<int dim>
bool Quadrature< dim >::operator== ( const Quadrature< dim > &  p) const

Test for equality of two quadratures.

Definition at line 263 of file quadrature.cc.

◆ initialize()

template<int dim>
void Quadrature< dim >::initialize ( const std::vector< Point< dim > > &  points,
const std::vector< double > &  weights 
)

Set the quadrature points and weights to the values provided in the arguments.

Definition at line 55 of file quadrature.cc.

◆ size()

template<int dim>
unsigned int Quadrature< dim >::size ( ) const

Number of quadrature points.

◆ point()

template<int dim>
const Point<dim>& Quadrature< dim >::point ( const unsigned int  i) const

Return the ith quadrature point.

◆ get_points()

template<int dim>
const std::vector<Point<dim> >& Quadrature< dim >::get_points ( ) const

Return a reference to the whole array of quadrature points.

◆ weight()

template<int dim>
double Quadrature< dim >::weight ( const unsigned int  i) const

Return the weight of the ith quadrature point.

◆ get_weights()

template<int dim>
const std::vector<double>& Quadrature< dim >::get_weights ( ) const

Return a reference to the whole array of weights.

◆ memory_consumption()

template<int dim>
std::size_t Quadrature< dim >::memory_consumption ( ) const

Determine an estimate for the memory consumption (in bytes) of this object.

Definition at line 280 of file quadrature.cc.

◆ serialize()

template<int dim>
template<class Archive >
void Quadrature< dim >::serialize ( Archive &  ar,
const unsigned int  version 
)

Write or read the data of this object to or from a stream for the purpose of serialization.

Member Data Documentation

◆ quadrature_points

template<int dim>
std::vector<Point<dim> > Quadrature< dim >::quadrature_points
protected

List of quadrature points. To be filled by the constructors of derived classes.

Definition at line 228 of file quadrature.h.

◆ weights

template<int dim>
std::vector<double> Quadrature< dim >::weights
protected

List of weights of the quadrature points. To be filled by the constructors of derived classes.

Definition at line 234 of file quadrature.h.


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