Intrepid
Public Member Functions | Private Attributes | List of all members
Intrepid::CubatureControlVolumeBoundary< Scalar, ArrayPoint, ArrayWeight > Class Template Reference

Defines cubature (integration) rules over Neumann boundaries for control volume method. More...

#include <Intrepid_CubatureControlVolumeBoundary.hpp>

Inheritance diagram for Intrepid::CubatureControlVolumeBoundary< Scalar, ArrayPoint, ArrayWeight >:
Intrepid::Cubature< Scalar, ArrayPoint, ArrayWeight >

Public Member Functions

 CubatureControlVolumeBoundary (const Teuchos::RCP< const shards::CellTopology > &cellTopology, int cellSide=0)
 
void getCubature (ArrayPoint &cubPoints, ArrayWeight &cubWeights) const
 Returns cubature points and weights Method for reference space cubature - throws an exception.
 
void getCubature (ArrayPoint &cubPoints, ArrayWeight &cubWeights, ArrayPoint &cellCoords) const
 Returns cubature points and weights (return arrays must be pre-sized/pre-allocated).
 
int getNumPoints () const
 Returns the number of cubature points.
 
int getDimension () const
 Returns dimension of integration domain.
 
void getAccuracy (std::vector< int > &accuracy) const
 Returns max. degree of polynomials that are integrated exactly. The return vector has size 1.
 
virtual void getCubature (ArrayPoint &cubPoints, ArrayWeight &cubWeights) const=0
 Returns cubature points and weights (return arrays must be pre-sized/pre-allocated).
 
virtual void getCubature (ArrayPoint &cubPoints, ArrayWeight &cubWeights, ArrayPoint &cellVertices) const=0
 Returns cubature points and weights on physical cells (return arrays must be pre-sized/pre-allocated).
 
virtual int getNumPoints () const=0
 Returns the number of cubature points.
 
virtual int getDimension () const=0
 Returns dimension of the integration domain.
 
virtual void getAccuracy (std::vector< int > &accuracy) const=0
 Returns algebraic accuracy (e.g. max. degree of polynomial that is integrated exactly). For tensor-product or sparse rules, algebraic accuracy for each coordinate direction is returned.
 

Private Attributes

Teuchos::RCP< const shards::CellTopology > primaryCellTopo_
 The topology of the primary cell side.
 
Teuchos::RCP< const shards::CellTopology > subCVCellTopo_
 The topology of the sub-control volume.
 
int degree_
 The degree of the polynomials that are integrated exactly.
 
int numPoints_
 The number of cubature points.
 
int cubDimension_
 Dimension of integration domain.
 
int sideIndex_
 Index of cell side.
 

Detailed Description

template<class Scalar, class ArrayPoint, class ArrayWeight>
class Intrepid::CubatureControlVolumeBoundary< Scalar, ArrayPoint, ArrayWeight >

Defines cubature (integration) rules over Neumann boundaries for control volume method.

Integration on Neumann boundaries for the control volume method requires integration points defined on primary cell sides. These points are not equivalent to control volume points on lower dimensional topologies and therefore require a separate class to define them.

Definition at line 69 of file Intrepid_CubatureControlVolumeBoundary.hpp.

Constructor & Destructor Documentation

◆ CubatureControlVolumeBoundary()

template<class Scalar , class ArrayPoint , class ArrayWeight >
Intrepid::CubatureControlVolumeBoundary< Scalar, ArrayPoint, ArrayWeight >::CubatureControlVolumeBoundary ( const Teuchos::RCP< const shards::CellTopology > &  cellTopology,
int  cellSide = 0 
)

brief Constructor.

Parameters
cellTopology[in] - The topology of the primary cell.
cellSide[in] - The index of the boundary side of the primary cell

Definition at line 55 of file Intrepid_CubatureControlVolumeBoundaryDef.hpp.

◆ ~CubatureControlVolumeBoundary()

template<class Scalar , class ArrayPoint , class ArrayWeight >
virtual Intrepid::CubatureControlVolumeBoundary< Scalar, ArrayPoint, ArrayWeight >::~CubatureControlVolumeBoundary ( )
inlinevirtual

Definition at line 113 of file Intrepid_CubatureControlVolumeBoundary.hpp.

Member Function Documentation

◆ getAccuracy()

template<class Scalar , class ArrayPoint , class ArrayWeight >
void Intrepid::CubatureControlVolumeBoundary< Scalar, ArrayPoint, ArrayWeight >::getAccuracy ( std::vector< int > &  accuracy) const
virtual

Returns max. degree of polynomials that are integrated exactly. The return vector has size 1.

Implements Intrepid::Cubature< Scalar, ArrayPoint, ArrayWeight >.

Definition at line 224 of file Intrepid_CubatureControlVolumeBoundaryDef.hpp.

◆ getCubature() [1/2]

template<class Scalar , class ArrayPoint , class ArrayWeight >
void Intrepid::CubatureControlVolumeBoundary< Scalar, ArrayPoint, ArrayWeight >::getCubature ( ArrayPoint &  cubPoints,
ArrayWeight &  cubWeights 
) const
virtual

Returns cubature points and weights Method for reference space cubature - throws an exception.

Parameters
cubPoints[out] - Array containing the cubature points.
cubWeights[out] - Array of corresponding cubature weights.

Implements Intrepid::Cubature< Scalar, ArrayPoint, ArrayWeight >.

Definition at line 79 of file Intrepid_CubatureControlVolumeBoundaryDef.hpp.

◆ getCubature() [2/2]

template<class Scalar , class ArrayPoint , class ArrayWeight >
void Intrepid::CubatureControlVolumeBoundary< Scalar, ArrayPoint, ArrayWeight >::getCubature ( ArrayPoint &  cubPoints,
ArrayWeight &  cubWeights,
ArrayPoint &  cellCoords 
) const
virtual

Returns cubature points and weights (return arrays must be pre-sized/pre-allocated).

Parameters
cubPoints[out] - Array containing the cubature points.
cubWeights[out] - Array of corresponding cubature weights.
cellCoords[in] - Array of cell coordinates

Implements Intrepid::Cubature< Scalar, ArrayPoint, ArrayWeight >.

Definition at line 87 of file Intrepid_CubatureControlVolumeBoundaryDef.hpp.

References Intrepid::CellTools< Scalar >::getSubCVCoords(), Intrepid::CellTools< Scalar >::mapToReferenceSubcell(), and Intrepid::CellTools< Scalar >::setJacobianDet().

◆ getDimension()

template<class Scalar , class ArrayPoint , class ArrayWeight >
int Intrepid::CubatureControlVolumeBoundary< Scalar, ArrayPoint, ArrayWeight >::getDimension
virtual

Returns dimension of integration domain.

Implements Intrepid::Cubature< Scalar, ArrayPoint, ArrayWeight >.

Definition at line 219 of file Intrepid_CubatureControlVolumeBoundaryDef.hpp.

◆ getNumPoints()

template<class Scalar , class ArrayPoint , class ArrayWeight >
int Intrepid::CubatureControlVolumeBoundary< Scalar, ArrayPoint, ArrayWeight >::getNumPoints
virtual

Returns the number of cubature points.

Implements Intrepid::Cubature< Scalar, ArrayPoint, ArrayWeight >.

Definition at line 214 of file Intrepid_CubatureControlVolumeBoundaryDef.hpp.

Member Data Documentation

◆ cubDimension_

template<class Scalar , class ArrayPoint , class ArrayWeight >
int Intrepid::CubatureControlVolumeBoundary< Scalar, ArrayPoint, ArrayWeight >::cubDimension_
private

Dimension of integration domain.

Definition at line 136 of file Intrepid_CubatureControlVolumeBoundary.hpp.

◆ degree_

template<class Scalar , class ArrayPoint , class ArrayWeight >
int Intrepid::CubatureControlVolumeBoundary< Scalar, ArrayPoint, ArrayWeight >::degree_
private

The degree of the polynomials that are integrated exactly.

Definition at line 128 of file Intrepid_CubatureControlVolumeBoundary.hpp.

◆ numPoints_

template<class Scalar , class ArrayPoint , class ArrayWeight >
int Intrepid::CubatureControlVolumeBoundary< Scalar, ArrayPoint, ArrayWeight >::numPoints_
private

The number of cubature points.

Definition at line 132 of file Intrepid_CubatureControlVolumeBoundary.hpp.

◆ primaryCellTopo_

template<class Scalar , class ArrayPoint , class ArrayWeight >
Teuchos::RCP<const shards::CellTopology> Intrepid::CubatureControlVolumeBoundary< Scalar, ArrayPoint, ArrayWeight >::primaryCellTopo_
private

The topology of the primary cell side.

Definition at line 120 of file Intrepid_CubatureControlVolumeBoundary.hpp.

◆ sideIndex_

template<class Scalar , class ArrayPoint , class ArrayWeight >
int Intrepid::CubatureControlVolumeBoundary< Scalar, ArrayPoint, ArrayWeight >::sideIndex_
private

Index of cell side.

Definition at line 140 of file Intrepid_CubatureControlVolumeBoundary.hpp.

◆ subCVCellTopo_

template<class Scalar , class ArrayPoint , class ArrayWeight >
Teuchos::RCP<const shards::CellTopology> Intrepid::CubatureControlVolumeBoundary< Scalar, ArrayPoint, ArrayWeight >::subCVCellTopo_
private

The topology of the sub-control volume.

Definition at line 124 of file Intrepid_CubatureControlVolumeBoundary.hpp.


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