Reference documentation for deal.II version GIT relicensing-136-gb80d0be4af 2024-03-18 08:20:02+00:00
\(\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
Classes | Public Types | Public Member Functions | Static Public Member Functions | List of all members
QProjector< dim > Class Template Reference

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

Classes

class  DataSetDescriptor
 

Public Types

using SubQuadrature = Quadrature< dim - 1 >
 

Public Member Functions

void project_to_face (const ReferenceCell &reference_cell, const Quadrature< 0 > &, 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)
 
Quadrature< 3 > project_to_oriented_face (const ReferenceCell &reference_cell, const Quadrature< 2 > &quadrature, const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation)
 
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 > &)
 
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< 3 > project_to_oriented_subface (const ReferenceCell &reference_cell, const Quadrature< 2 > &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< 3 > ref_case)
 
Quadrature< 1 > project_to_all_faces (const ReferenceCell &reference_cell, const hp::QCollection< 0 > &quadrature)
 
Quadrature< 2 > project_to_all_faces (const ReferenceCell &reference_cell, const hp::QCollection< 1 > &quadrature)
 
Quadrature< 3 > project_to_all_faces (const ReferenceCell &reference_cell, const hp::QCollection< 2 > &quadrature)
 
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 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_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)
 

Detailed Description

template<int dim>
class QProjector< dim >

This class is a helper class to facilitate the usage of quadrature formulae on faces or subfaces of cells. It computes the locations of quadrature points on the unit cell from a quadrature object for a manifold of one dimension less than that of the cell and the number of the face. For example, giving the Simpson rule in one dimension and using the project_to_face() function with face number 1, the returned points will be (1,0), (1,0.5) and (1,1). Note that faces have an orientation, so when projecting to face 3, you will get (0,0), (0,0.5) and (0,1), which is in clockwise sense, while for face 1 the points were in counterclockwise sense.

For the projection to subfaces (i.e. to the children of a face of the unit cell), the same applies as above. Note the order in which the children of a face are numbered, which in two dimensions coincides with the orientation of the face.

The second set of functions generates a quadrature formula by projecting a given quadrature rule on all faces and subfaces. This is used in the FEFaceValues and FESubfaceValues classes. Since we now have the quadrature points of all faces and subfaces in one array, we need to have a way to find the starting index of the points and weights corresponding to one face or subface within this array. This is done through the DataSetDescriptor member class.

The different functions are grouped into a common class to avoid putting them into global namespace. However, since they have no local data, all functions are declared static and can be called without creating an object of this class.

For the 3d case, you should note that the orientation of faces is even more intricate than for two dimensions. Quadrature formulae are projected upon the faces in their standard orientation, not to the inside or outside of the hexahedron. To make things more complicated, in 3d we allow faces in two orientations (which can be identified using cell->face_orientation(face)), so we have to project quadrature formula onto faces and subfaces in two orientations. (Refer to the documentation of the Triangulation class for a description of the orientation of the different faces, as well as to the glossary entry on face orientation for more information on this.) The DataSetDescriptor member class is used to identify where each dataset starts.

Definition at line 81 of file qprojector.h.

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 88 of file qprojector.h.

Member Function Documentation

◆ project_to_face() [1/5]

template<int dim>
static 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.

◆ project_to_face() [2/5]

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.

Definition at line 1721 of file qprojector.cc.

◆ project_to_oriented_face() [1/2]

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.

Definition at line 442 of file qprojector.cc.

◆ project_to_subface() [1/5]

template<int dim>
static 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.

◆ project_to_subface() [2/5]

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.

Definition at line 1738 of file qprojector.cc.

◆ project_to_oriented_subface() [1/2]

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.

Definition at line 664 of file qprojector.cc.

◆ project_to_all_faces() [1/5]

template<int dim>
static 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
In 3d, this function produces eight sets of quadrature points for each face, in order to cope possibly different orientations of the mesh.

◆ project_to_all_faces() [2/5]

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 409 of file qprojector.h.

◆ project_to_all_subfaces() [1/4]

template<int dim>
static 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.

Definition at line 1226 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.

Definition at line 1258 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 1290 of file qprojector.cc.

◆ project_to_face() [3/5]

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

Definition at line 341 of file qprojector.cc.

◆ project_to_face() [4/5]

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 360 of file qprojector.cc.

◆ project_to_face() [5/5]

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 422 of file qprojector.cc.

◆ project_to_oriented_face() [2/2]

Quadrature< 3 > QProjector< 3 >::project_to_oriented_face ( const ReferenceCell reference_cell,
const Quadrature< 2 > &  quadrature,
const unsigned int  face_no,
const bool  face_orientation,
const bool  face_flip,
const bool  face_rotation 
)

Definition at line 456 of file qprojector.cc.

◆ project_to_subface() [3/5]

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 475 of file qprojector.cc.

◆ project_to_subface() [4/5]

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 > &   
)

Definition at line 496 of file qprojector.cc.

◆ project_to_subface() [5/5]

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 640 of file qprojector.cc.

◆ project_to_oriented_subface() [2/2]

Quadrature< 3 > QProjector< 3 >::project_to_oriented_subface ( const ReferenceCell reference_cell,
const Quadrature< 2 > &  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< 3 >  ref_case 
)

Definition at line 686 of file qprojector.cc.

◆ project_to_all_faces() [3/5]

Quadrature< 1 > QProjector< 1 >::project_to_all_faces ( const ReferenceCell reference_cell,
const hp::QCollection< 0 > &  quadrature 
)

Definition at line 718 of file qprojector.cc.

◆ project_to_all_faces() [4/5]

Quadrature< 2 > QProjector< 2 >::project_to_all_faces ( const ReferenceCell reference_cell,
const hp::QCollection< 1 > &  quadrature 
)

Definition at line 765 of file qprojector.cc.

◆ project_to_all_faces() [5/5]

Quadrature< 3 > QProjector< 3 >::project_to_all_faces ( const ReferenceCell reference_cell,
const hp::QCollection< 2 > &  quadrature 
)

Definition at line 890 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 1062 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 1109 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 1160 of file qprojector.cc.


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