deal.II version GIT relicensing-2165-gc91f007519 2024-11-20 01:40:00+00:00
|
#include <deal.II/base/quadrature.h>
Public Types | |
using | SubQuadrature = Quadrature< dim==0 ? 0 :dim - 1 > |
Public Member Functions | |
Quadrature () | |
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 (std::vector< Point< dim > > &&points, std::vector< double > &&weights) | |
Quadrature (const std::vector< Point< dim > > &points) | |
Quadrature (const Point< dim > &point) | |
virtual | ~Quadrature () override=default |
Quadrature & | operator= (const Quadrature< dim > &) |
Quadrature & | operator= (Quadrature< dim > &&)=default |
bool | operator== (const Quadrature< dim > &p) const |
void | initialize (const ArrayView< const Point< dim > > &points, const ArrayView< const double > &weights={}) |
unsigned int | size () const |
bool | empty () 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 |
EnableObserverPointer functionality | |
Classes derived from EnableObserverPointer provide a facility to subscribe to this object. This is mostly used by the ObserverPointer 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 |
Static Public Member Functions | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Protected Member Functions | |
Quadrature (const unsigned int n_quadrature_points) | |
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 Types | |
using | map_value_type = decltype(counter_map)::value_type |
using | map_iterator = decltype(counter_map)::iterator |
Private Member Functions | |
void | check_no_subscribers () const noexcept |
Private Attributes | |
std::atomic< unsigned int > | counter |
std::map< std::string, unsigned int > | counter_map |
std::vector< std::atomic< bool > * > | validity_pointers |
const std::type_info * | object_info |
Static Private Attributes | |
static std::mutex | mutex |
Base class for quadrature formulae in arbitrary dimensions. Quadrature is a means to approximate an integral by evaluating the integrand at specific points \(\mathbf x_q\) and summing the point values with specific weights \(w_q\); that is, quadrature computes
\begin{align*} \int_K f(\mathbf x) \; dx \approx \sum_{q=0,\ldots,Q-1} f(\mathbf x_q) w_q. \end{align*}
This class stores quadrature points \(\mathbf x_q\) and weights \(w_q\) for concrete "quadrature formulas" when \(K\) (the domain we integrate over) is a reference cell. That is, points and weights are expressed in the coordinate system of a reference cell (see the ReferenceCell class) and as such serves to represent quadrature points and weights on the unit line segment \([0,1]\) in 1d, on the unit square or unit triangle in 2d, as well as the unit tetrahedron, cube, pyramid, and wedge reference cells in 3d. Integration over concrete cells is done by coordinate transformation to the reference cell represented by the current class.
There are a number of derived classes, denoting concrete integration formulae. Their names are prefixed by Q
. Refer to the list of derived classes for more details.
At least for quadrilaterals and hexahedra (or, more precisely, since we work on reference cells: for the unit square and the unit cube), quadrature formulas are typically tensor products of one-dimensional formulas (see also 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.
For each quadrature formula we denote by m
, the maximal degree of polynomials integrated exactly on the reference cell the quadrature formula corresponds to. 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+1)
st order.
At least for hypercube reference cells (i.e., squares and cubes), 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.
Quadrature objects are used in a number of places within deal.II where integration is performed, most notably via the FEValues and related classes. Some of these classes are also used in contexts where no integrals are involved, but where functions need to be evaluated at specific points, for example to evaluate the solution at individual points or to create graphical output. Examples are the implementation of VectorTools::point_value() and the DataOut and related classes (in particular in connection with the DataPostprocessor class). In such contexts, one often creates specific "Quadrature" objects in which the "quadrature points" are simply the points (in the coordinate system of the reference cell) at which one wants to evaluate the solution. In these kinds of cases, the weights stored by the current class are not used and the name "quadrature object" is interpreted as "list of evaluation points".
Definition at line 122 of file quadrature.h.
using Quadrature< dim >::SubQuadrature = Quadrature<dim == 0 ? 0 : dim - 1> |
Define an alias for a quadrature that acts on an object of one dimension less. For cells, this would then be a face quadrature. A sub quadrature of a 0-dimensional quadrature is defined as still being 0-dimensional.
Definition at line 130 of file quadrature.h.
|
privateinherited |
The data type used in counter_map.
Definition at line 238 of file enable_observer_pointer.h.
|
privateinherited |
The iterator type used in counter_map.
Definition at line 243 of file enable_observer_pointer.h.
Quadrature< dim >::Quadrature | ( | ) |
Default constructor.
Definition at line 46 of file quadrature.cc.
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. This constructor assumes (and tests) that constant functions are integrated exactly, i.e. the sum of the quadrature weights is one.
SubQuadrature<dim>::type
expands to Quadrature<dim-1>
.
Definition at line 169 of file quadrature.cc.
|
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.
This constructor does not require that constant functions are integrated exactly. Therefore, it is appropriate if the one-dimensional formula is defined with respect to a weighting function.
If dim == 0, the resulting quadrature formula will be a single Point<0> having unit weight.
Definition at line 267 of file quadrature.cc.
Quadrature< dim >::Quadrature | ( | const Quadrature< dim > & | q | ) |
Copy constructor.
Definition at line 305 of file quadrature.cc.
|
defaultnoexcept |
Move constructor. Construct a new quadrature object by transferring the internal data of another quadrature object.
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 87 of file quadrature.cc.
Quadrature< dim >::Quadrature | ( | std::vector< Point< dim > > && | points, |
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, moving the points and weights into the present object.
Definition at line 100 of file quadrature.cc.
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 113 of file quadrature.cc.
Quadrature< dim >::Quadrature | ( | const Point< dim > & | point | ) |
Constructor for a one-point quadrature. Sets the weight of this point to one.
Definition at line 125 of file quadrature.cc.
|
overridevirtualdefault |
Virtual destructor.
|
explicitprotected |
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. Nonetheless, it is easy to accidentally write
where QGauss was meant. As a consequence, this constructor is protected
and so only available to derived classes initializing their base class.
Definition at line 53 of file quadrature.cc.
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 320 of file quadrature.cc.
|
default |
Move assignment operator. Moves all data from another quadrature object to this object.
bool Quadrature< dim >::operator== | ( | const Quadrature< dim > & | p | ) | const |
Test for equality of two quadratures.
Definition at line 340 of file quadrature.cc.
void Quadrature< dim >::initialize | ( | const ArrayView< const Point< dim > > & | points, |
const ArrayView< const double > & | weights = {} |
||
) |
Set the quadrature points and weights to the values provided in the arguments. The weights array is allowed to be empty, in which case the weights are 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 63 of file quadrature.cc.
unsigned int Quadrature< dim >::size | ( | ) | const |
Number of quadrature points.
bool Quadrature< dim >::empty | ( | ) | const |
Return if quadrature is empty.
const Point< dim > & Quadrature< dim >::point | ( | const unsigned int | i | ) | const |
Return the i
th quadrature point.
const std::vector< Point< dim > > & Quadrature< dim >::get_points | ( | ) | const |
Return a reference to the whole array of quadrature points.
double Quadrature< dim >::weight | ( | const unsigned int | i | ) | const |
Return the weight of the i
th quadrature point.
const std::vector< double > & Quadrature< dim >::get_weights | ( | ) | const |
Return a reference to the whole array of weights.
std::size_t Quadrature< dim >::memory_consumption | ( | ) | const |
Determine an estimate for the memory consumption (in bytes) of this object.
Definition at line 349 of file quadrature.cc.
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 using the BOOST serialization library.
bool Quadrature< dim >::is_tensor_product | ( | ) | const |
This function returns true if the quadrature object is a tensor product of one-dimensional formulas and the quadrature points are sorted lexicographically.
std::conditional_t< dim==1, std::array< Quadrature< 1 >, dim >, const std::array< Quadrature< 1 >, dim > & > Quadrature< dim >::get_tensor_basis | ( | ) | const |
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
).
Definition at line 361 of file quadrature.cc.
|
inherited |
Subscribes a user of the object by storing the pointer validity
. The subscriber may be identified by text supplied as identifier
.
Definition at line 131 of file enable_observer_pointer.cc.
|
inherited |
Unsubscribes a user from the object.
identifier
and the validity
pointer must be the same as the one supplied to subscribe(). Definition at line 151 of file enable_observer_pointer.cc.
|
inlineinherited |
Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.
Definition at line 322 of file enable_observer_pointer.h.
|
inlineinherited |
List the subscribers to the input stream
.
Definition at line 339 of file enable_observer_pointer.h.
|
inherited |
List the subscribers to deallog
.
Definition at line 199 of file enable_observer_pointer.cc.
|
privatenoexceptinherited |
Check that there are no objects subscribing to this object. If this check passes then it is safe to destroy the current object. It this check fails then this function will either abort or print an error message to deallog (by using the AssertNothrow mechanism), but will not throw an exception.
Definition at line 53 of file enable_observer_pointer.cc.
|
protected |
List of quadrature points. To be filled by the constructors of derived classes.
Definition at line 353 of file quadrature.h.
|
protected |
List of weights of the quadrature points. To be filled by the constructors of derived classes.
Definition at line 359 of file quadrature.h.
|
protected |
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 368 of file quadrature.h.
|
protected |
Stores the one-dimensional tensor basis objects in case this object can be represented by a tensor product.
Definition at line 374 of file quadrature.h.
|
mutableprivateinherited |
Store the number of objects which subscribed to this object. Initially, this number is zero, and upon destruction it shall be zero again (i.e. all objects which subscribed should have unsubscribed again).
The creator (and owner) of an object is counted in the map below if HE manages to supply identification.
We use the mutable
keyword in order to allow subscription to constant objects also.
This counter may be read from and written to concurrently in multithreaded code: hence we use the std::atomic
class template.
Definition at line 227 of file enable_observer_pointer.h.
|
mutableprivateinherited |
In this map, we count subscriptions for each different identification string supplied to subscribe().
Definition at line 233 of file enable_observer_pointer.h.
|
mutableprivateinherited |
In this vector, we store pointers to the validity bool in the ObserverPointer objects that subscribe to this class.
Definition at line 249 of file enable_observer_pointer.h.
|
mutableprivateinherited |
Pointer to the typeinfo object of this object, from which we can later deduce the class name. Since this information on the derived class is neither available in the destructor, nor in the constructor, we obtain it in between and store it here.
Definition at line 257 of file enable_observer_pointer.h.
|
staticprivateinherited |
A mutex used to ensure data consistency when accessing the mutable
members of this class. This lock is used in the subscribe() and unsubscribe() functions, as well as in list_subscribers()
.
Definition at line 280 of file enable_observer_pointer.h.