deal.II version 9.7.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
QProjector< dim > Class Template Reference

Class which transforms dim - 1-dimensional quadrature rules to dim-dimensional face quadratures. More...

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

Detailed Description

template<int dim>
class QProjector< dim >

Class which transforms dim - 1-dimensional quadrature rules to dim-dimensional face quadratures.

The majority of the finite element infrastructure, such as FE_Q and FE_SimplexP, uses polynomials defined on a reference cell: for example, FE_Q<3> defines a polynomial space whose domain is the unit hexahedron (i.e., ReferenceCells::Hexahedron). Hence, computing quadratures using shape functions on a face of a reference cell requires converting a lower-dimensional Quadrature into one defined on the boundary of the higher-dimensional object, e.g., converting a Quadrature defined on a quadrilateral into one defined on one face of a hexahedron.

QProjector computes the locations of quadrature points on faces or subfaces of reference cells from a Quadrature of one dimension less than that of the cell, face (and possibly also subface) number, and orientation. For example, calling QProjector::project_to_face() with QSimpson<1>, face number 1, and numbers::default_geometric_orientation returns a Quadrature with points (1,0), (1,0.5), and (1,1) with weights equal to the 1d case. Similarly, if we instead use face 3 and numbers::reverse_line_orientation we obtain points (1,1), (0.5,1) and (0,1). Projection to subfaces works in the same way.

In practice, computing face integrals (e.g., via FEFaceValues or FESubfaceValues) requires quadrature rules for all possible permutations of face number, subface number, and combined orientation. This class provides several functions for doing just that, such as QProjector::project_to_all_faces(). Furthermore, the DataSetDescriptor class implements indexing for converting between face number, subface number, and combined orientation to the index of the associated Quadrature rule.

Definition at line 68 of file qprojector.h.

Classes

class  DataSetDescriptor
 Class storing the offset index into a Quadrature rule created by project_to_all_faces() or project_to_all_subfaces(). More...

Public Types

using SubQuadrature = Quadrature<dim - 1>

Public Member Functions

void project_to_face (const ReferenceCell &reference_cell, const Quadrature< 0 > &quadrature, const unsigned int face_no, std::vector< Point< 1 > > &q_points)
void project_to_face (const ReferenceCell &reference_cell, const Quadrature< 1 > &quadrature, const unsigned int face_no, std::vector< Point< 2 > > &q_points)
void project_to_face (const ReferenceCell &reference_cell, const Quadrature< 2 > &quadrature, const unsigned int face_no, std::vector< Point< 3 > > &q_points)
void project_to_subface (const ReferenceCell &reference_cell, const Quadrature< 0 > &, const unsigned int face_no, const unsigned int, std::vector< Point< 1 > > &q_points, const RefinementCase< 0 > &)
void project_to_subface (const ReferenceCell &reference_cell, const Quadrature< 1 > &quadrature, const unsigned int face_no, const unsigned int subface_no, std::vector< Point< 2 > > &q_points, const RefinementCase< 1 > &ref_case)
void project_to_subface (const ReferenceCell &reference_cell, const Quadrature< 2 > &quadrature, const unsigned int face_no, const unsigned int subface_no, std::vector< Point< 3 > > &q_points, const RefinementCase< 2 > &ref_case)
Quadrature< 1 > project_to_all_subfaces (const ReferenceCell &reference_cell, const Quadrature< 0 > &quadrature)
Quadrature< 2 > project_to_all_subfaces (const ReferenceCell &reference_cell, const SubQuadrature &quadrature)
Quadrature< 3 > project_to_all_subfaces (const ReferenceCell &reference_cell, const SubQuadrature &quadrature)

Static Public Member Functions

static void project_to_face (const ReferenceCell &reference_cell, const SubQuadrature &quadrature, const unsigned int face_no, std::vector< Point< dim > > &q_points)
static Quadrature< dim > project_to_face (const ReferenceCell &reference_cell, const SubQuadrature &quadrature, const unsigned int face_no)
static Quadrature< dim > project_to_oriented_face (const ReferenceCell &reference_cell, const SubQuadrature &quadrature, const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation)
static Quadrature< dim > project_to_face (const ReferenceCell &reference_cell, const SubQuadrature &quadrature, const unsigned int face_no, const types::geometric_orientation combined_orientation)
static void project_to_subface (const ReferenceCell &reference_cell, const SubQuadrature &quadrature, const unsigned int face_no, const unsigned int subface_no, std::vector< Point< dim > > &q_points, const RefinementCase< dim - 1 > &ref_case=RefinementCase< dim - 1 >::isotropic_refinement)
static Quadrature< dim > project_to_subface (const ReferenceCell &reference_cell, const SubQuadrature &quadrature, const unsigned int face_no, const unsigned int subface_no, const RefinementCase< dim - 1 > &ref_case=RefinementCase< dim - 1 >::isotropic_refinement)
static Quadrature< dim > project_to_oriented_subface (const ReferenceCell &reference_cell, const SubQuadrature &quadrature, const unsigned int face_no, const unsigned int subface_no, const bool face_orientation, const bool face_flip, const bool face_rotation, const internal::SubfaceCase< dim > ref_case)
static Quadrature< dim > project_to_subface (const ReferenceCell &reference_cell, const SubQuadrature &quadrature, const unsigned int face_no, const unsigned int subface_no, const types::geometric_orientation combined_orientation, const RefinementCase< dim - 1 > &ref_case)
static Quadrature< dim > project_to_all_faces (const ReferenceCell &reference_cell, const hp::QCollection< dim - 1 > &quadrature)
static Quadrature< dim > project_to_all_faces (const ReferenceCell &reference_cell, const Quadrature< dim - 1 > &quadrature)
static Quadrature< dim > project_to_all_subfaces (const ReferenceCell &reference_cell, const SubQuadrature &quadrature)
static Quadrature< dim > project_to_child (const ReferenceCell &reference_cell, const Quadrature< dim > &quadrature, const unsigned int child_no)
static Quadrature< dim > project_to_all_children (const ReferenceCell &reference_cell, const Quadrature< dim > &quadrature)
static Quadrature< dim > project_to_line (const ReferenceCell &reference_cell, const Quadrature< 1 > &quadrature, const Point< dim > &p1, const Point< dim > &p2)

Member Typedef Documentation

◆ SubQuadrature

template<int dim>
using QProjector< dim >::SubQuadrature = Quadrature<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.

Definition at line 75 of file qprojector.h.

Member Function Documentation

◆ project_to_face() [1/6]

template<int dim>
void QProjector< dim >::project_to_face ( const ReferenceCell & reference_cell,
const SubQuadrature & quadrature,
const unsigned int face_no,
std::vector< Point< dim > > & q_points )
static

Compute the quadrature points on the cell if the given quadrature formula is used on face face_no. For further details, see the general doc for this class.

Deprecated
Use the version of this function which takes a combined_orientation argument instead.

◆ project_to_face() [2/6]

template<int dim>
Quadrature< dim > QProjector< dim >::project_to_face ( const ReferenceCell & reference_cell,
const SubQuadrature & quadrature,
const unsigned int face_no )
static

Compute the cell quadrature formula corresponding to using quadrature on face face_no. For further details, see the general doc for this class.

Deprecated
Use the version of this function which takes a combined_orientation argument instead.

Definition at line 1104 of file qprojector.cc.

◆ project_to_oriented_face()

template<int dim>
Quadrature< dim > QProjector< dim >::project_to_oriented_face ( const ReferenceCell & reference_cell,
const SubQuadrature & quadrature,
const unsigned int face_no,
const bool face_orientation,
const bool face_flip,
const bool face_rotation )
static

Compute the cell quadrature formula corresponding to using quadrature on face face_no taking into account the orientation of the face. For further details, see the general doc for this class.

Deprecated
Use the version of project_to_face() which takes a combined_orientation argument instead.

Definition at line 319 of file qprojector.cc.

◆ project_to_face() [3/6]

template<int dim>
Quadrature< dim > QProjector< dim >::project_to_face ( const ReferenceCell & reference_cell,
const SubQuadrature & quadrature,
const unsigned int face_no,
const types::geometric_orientation combined_orientation )
static

Compute the cell quadrature formula corresponding to using quadrature on face face_no. For further details, see the general doc for this class.

Definition at line 339 of file qprojector.cc.

◆ project_to_subface() [1/6]

template<int dim>
void QProjector< dim >::project_to_subface ( const ReferenceCell & reference_cell,
const SubQuadrature & quadrature,
const unsigned int face_no,
const unsigned int subface_no,
std::vector< Point< dim > > & q_points,
const RefinementCase< dim - 1 > & ref_case = RefinementCase< dim - 1 >::isotropic_refinement )
static

Compute the quadrature points on the cell if the given quadrature formula is used on face face_no, subface number subface_no corresponding to RefineCase::Type ref_case. The last argument is only used in 3d.

Note
Only the points are transformed. The quadrature weights are the same as those of the original rule.
Deprecated
Use the version of project_to_subface() which takes a combined_orientation argument instead.

◆ project_to_subface() [2/6]

template<int dim>
Quadrature< dim > QProjector< dim >::project_to_subface ( const ReferenceCell & reference_cell,
const SubQuadrature & quadrature,
const unsigned int face_no,
const unsigned int subface_no,
const RefinementCase< dim - 1 > & ref_case = RefinementCase<dim - 1>::isotropic_refinement )
static

Compute the cell quadrature formula corresponding to using quadrature on subface subface_no of face face_no with RefinementCase<dim-1> ref_case. The last argument is only used in 3d.

Note
Only the points are transformed. The quadrature weights are the same as those of the original rule.
This function is deprecated since it makes an implicit assumption that the cell is a line (1D), a quad (2d), or a hex (3d). Use the other version of this function that takes the reference cell type instead.
Deprecated
Use the version of project_to_subface() which takes a combined_orientation argument instead.

Definition at line 1121 of file qprojector.cc.

◆ project_to_oriented_subface()

template<int dim>
Quadrature< dim > QProjector< dim >::project_to_oriented_subface ( const ReferenceCell & reference_cell,
const SubQuadrature & quadrature,
const unsigned int face_no,
const unsigned int subface_no,
const bool face_orientation,
const bool face_flip,
const bool face_rotation,
const internal::SubfaceCase< dim > ref_case )
static

Compute the cell quadrature formula corresponding to using quadrature on subface subface_no of face face_no with SubfaceCase<dim> ref_case. The last argument is only used in 3d.

Note
Only the points are transformed. The quadrature weights are the same as those of the original rule.
Deprecated
Use the version of project_to_subface() which takes a combined_orientation argument instead.

Definition at line 442 of file qprojector.cc.

◆ project_to_subface() [3/6]

template<int dim>
Quadrature< dim > QProjector< dim >::project_to_subface ( const ReferenceCell & reference_cell,
const SubQuadrature & quadrature,
const unsigned int face_no,
const unsigned int subface_no,
const types::geometric_orientation combined_orientation,
const RefinementCase< dim - 1 > & ref_case )
static

Compute the cell quadrature formula corresponding to using quadrature on subface subface_no of face face_no with RefinementCase<dim-1> ref_case. The last argument is only used in 3d.

Note
Only the points are transformed. The quadrature weights are the same as those of the original rule.

Definition at line 464 of file qprojector.cc.

◆ project_to_all_faces() [1/2]

template<int dim>
Quadrature< dim > QProjector< dim >::project_to_all_faces ( const ReferenceCell & reference_cell,
const hp::QCollection< dim - 1 > & quadrature )
static

Take a collection of face quadrature formulas and generate a cell quadrature formula from it where the quadrature points of the given argument are projected on all faces.

The weights of the new rule are replications of the original weights. Thus, the sum of the weights is not one, but the number of faces, which is the surface of the reference cell.

This in particular allows us to extract a subset of points corresponding to a single face and use it as a quadrature on this face, as is done in FEFaceValues.

Note
This function creates ReferenceCell::n_face_orientations() sets of quadrature points for each face which are indexed (by combined orientation and face number) by a DataSetDescriptor.

Definition at line 585 of file qprojector.cc.

◆ project_to_all_faces() [2/2]

template<int dim>
Quadrature< dim > QProjector< dim >::project_to_all_faces ( const ReferenceCell & reference_cell,
const Quadrature< dim - 1 > & quadrature )
inlinestatic

Like the above function, applying the same face quadrature formula on all faces.

Definition at line 532 of file qprojector.h.

◆ project_to_all_subfaces() [1/4]

template<int dim>
Quadrature< dim > QProjector< dim >::project_to_all_subfaces ( const ReferenceCell & reference_cell,
const SubQuadrature & quadrature )
static

Take a face quadrature formula and generate a cell quadrature formula from it where the quadrature points of the given argument are projected on all subfaces.

Like in project_to_all_faces(), the weights of the new rule sum up to the number of faces (not subfaces), which is the surface of the reference cell.

This in particular allows us to extract a subset of points corresponding to a single subface and use it as a quadrature on this face, as is done in FESubfaceValues.

◆ project_to_child()

template<int dim>
Quadrature< dim > QProjector< dim >::project_to_child ( const ReferenceCell & reference_cell,
const Quadrature< dim > & quadrature,
const unsigned int child_no )
static

Project a given quadrature formula to a child of a cell. You may want to use this function in case you want to extend an integral only over the area which a potential child would occupy. The child numbering is the same as the children would be numbered upon refinement of the cell.

As integration using this quadrature formula now only extends over a fraction of the cell, the weights of the resulting object are divided by GeometryInfo<dim>::children_per_cell.

Warning
This function is only implemented for hypercube elements.

Definition at line 782 of file qprojector.cc.

◆ project_to_all_children()

template<int dim>
Quadrature< dim > QProjector< dim >::project_to_all_children ( const ReferenceCell & reference_cell,
const Quadrature< dim > & quadrature )
static

Project a quadrature rule to all children of a cell. Similarly to project_to_all_subfaces(), this function replicates the formula generated by project_to_child() for all children, such that the weights sum up to one, the volume of the total cell again.

The child numbering is the same as the children would be numbered upon refinement of the cell.

Warning
This function is only implemented for hypercube elements.

Definition at line 814 of file qprojector.cc.

◆ project_to_line()

template<int dim>
Quadrature< dim > QProjector< dim >::project_to_line ( const ReferenceCell & reference_cell,
const Quadrature< 1 > & quadrature,
const Point< dim > & p1,
const Point< dim > & p2 )
static

Project the one dimensional rule quadrature to the straight line connecting the points p1 and p2.

Definition at line 846 of file qprojector.cc.

◆ project_to_face() [4/6]

void QProjector< 1 >::project_to_face ( const ReferenceCell & reference_cell,
const Quadrature< 0 > & quadrature,
const unsigned int face_no,
std::vector< Point< 1 > > & q_points )

Definition at line 262 of file qprojector.cc.

◆ project_to_face() [5/6]

void QProjector< 2 >::project_to_face ( const ReferenceCell & reference_cell,
const Quadrature< 1 > & quadrature,
const unsigned int face_no,
std::vector< Point< 2 > > & q_points )

Definition at line 280 of file qprojector.cc.

◆ project_to_face() [6/6]

void QProjector< 3 >::project_to_face ( const ReferenceCell & reference_cell,
const Quadrature< 2 > & quadrature,
const unsigned int face_no,
std::vector< Point< 3 > > & q_points )

Definition at line 298 of file qprojector.cc.

◆ project_to_subface() [4/6]

void QProjector< 1 >::project_to_subface ( const ReferenceCell & reference_cell,
const Quadrature< 0 > & ,
const unsigned int face_no,
const unsigned int ,
std::vector< Point< 1 > > & q_points,
const RefinementCase< 0 > &  )

Definition at line 375 of file qprojector.cc.

◆ project_to_subface() [5/6]

void QProjector< 2 >::project_to_subface ( const ReferenceCell & reference_cell,
const Quadrature< 1 > & quadrature,
const unsigned int face_no,
const unsigned int subface_no,
std::vector< Point< 2 > > & q_points,
const RefinementCase< 1 > & ref_case )

Definition at line 396 of file qprojector.cc.

◆ project_to_subface() [6/6]

void QProjector< 3 >::project_to_subface ( const ReferenceCell & reference_cell,
const Quadrature< 2 > & quadrature,
const unsigned int face_no,
const unsigned int subface_no,
std::vector< Point< 3 > > & q_points,
const RefinementCase< 2 > & ref_case )

Definition at line 418 of file qprojector.cc.

◆ project_to_all_subfaces() [2/4]

Quadrature< 1 > QProjector< 1 >::project_to_all_subfaces ( const ReferenceCell & reference_cell,
const Quadrature< 0 > & quadrature )

Definition at line 621 of file qprojector.cc.

◆ project_to_all_subfaces() [3/4]

Quadrature< 2 > QProjector< 2 >::project_to_all_subfaces ( const ReferenceCell & reference_cell,
const SubQuadrature & quadrature )

Definition at line 668 of file qprojector.cc.

◆ project_to_all_subfaces() [4/4]

Quadrature< 3 > QProjector< 3 >::project_to_all_subfaces ( const ReferenceCell & reference_cell,
const SubQuadrature & quadrature )

Definition at line 716 of file qprojector.cc.


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