Intrepid2
Public Types | Public Member Functions | Public Attributes | Static Public Attributes | List of all members
Intrepid2::Experimental::ProjectionStruct< SpT, ValueType > Class Template Reference

An helper class to compute the evaluation points and weights needed for performing projections. More...

#include <Intrepid2_ProjectionStruct.hpp>

Public Types

enum  EvalPointsType { BASIS, TARGET }
 
typedef Kokkos::pair
< ordinal_type, ordinal_type > 
range_type
 
typedef Kokkos::Impl::is_space
< SpT >
::host_mirror_space::execution_space 
host_space_type
 
typedef Kokkos::DynRankView
< ValueType, host_space_type > 
view_type
 
typedef Kokkos::View
< range_type
**, host_space_type > 
range_tag
 
typedef std::array< std::array
< view_type, maxSubCellsCount >
, numberSubCellDims > 
view_tag
 
typedef Kokkos::View< unsigned
**, host_space_type > 
key_tag
 

Public Member Functions

ordinal_type getNumBasisEvalPoints ()
 Returns number of basis evaluation points.
 
ordinal_type getNumBasisDerivEvalPoints ()
 Returns number of evaluation points for basis derivatives.
 
ordinal_type getNumTargetEvalPoints ()
 Returns number of points where to evaluate the target function.
 
ordinal_type getNumTargetDerivEvalPoints ()
 Returns number of points where to evaluate the derivatives of the target function.
 
ordinal_type getMaxNumDerivPoints (const EvalPointsType type) const
 Returns the maximum number of derivative evaluation points across all the subcells. More...
 
ordinal_type getMaxNumEvalPoints (const EvalPointsType type) const
 Returns the maximum number of evaluation points across all the subcells. More...
 
view_type getBasisEvalPoints (const ordinal_type subCellDim, const ordinal_type subCellId)
 Returns the basis evaluation points on a subcell. More...
 
view_type getBasisDerivEvalPoints (const ordinal_type subCellDim, const ordinal_type subCellId)
 Returns the evaluation points for basis derivatives on a subcell. More...
 
view_type getTargetEvalPoints (const ordinal_type subCellDim, const ordinal_type subCellId)
 Returns the points where to evaluate the target function on a subcell. More...
 
view_type getTargetDerivEvalPoints (const ordinal_type subCellDim, const ordinal_type subCellId)
 Returns the points where to evaluate the derivatives of the target function on a subcell. More...
 
view_type getEvalPoints (const ordinal_type subCellDim, const ordinal_type subCellId, EvalPointsType type) const
 Returns the basis/target evaluation points on a subcell. More...
 
view_type getDerivEvalPoints (const ordinal_type subCellDim, const ordinal_type subCellId, EvalPointsType type) const
 Returns the evaluation points for basis/target derivatives on a subcell. More...
 
view_type getBasisEvalWeights (const ordinal_type subCellDim, const ordinal_type subCellId)
 Returns the basis evaluation weights on a subcell. More...
 
view_type getBasisDerivEvalWeights (const ordinal_type subCellDim, const ordinal_type subCellId)
 Returns the basis derivatives evaluation weights on a subcell. More...
 
view_type getTargetEvalWeights (const ordinal_type subCellDim, const ordinal_type subCellId)
 Returns the function evaluation weights on a subcell. More...
 
view_type getTargetDerivEvalWeights (const ordinal_type subCellDim, const ordinal_type subCellId)
 Returns the function derivatives evaluation weights on a subcell. More...
 
const range_tag getBasisPointsRange () const
 Returns the range tag of the basis evaluation points subcells. More...
 
const range_tag getPointsRange (const EvalPointsType type) const
 Returns the range tag of the basis/target evaluation points in subcells. More...
 
const range_tag getBasisDerivPointsRange () const
 Returns the range tag of the derivative evaluation points on subcell. More...
 
const range_tag getDerivPointsRange (const EvalPointsType type) const
 Returns the range tag of the basis/target derivative evaluation points on subcells. More...
 
const range_tag getTargetPointsRange () const
 Returns the range of the target function evaluation points on subcells. More...
 
const range_tag getTargetDerivPointsRange () const
 Returns the range tag of the target function derivative evaluation points on subcells. More...
 
const key_tag getTopologyKey () const
 Returns the key tag for subcells. More...
 
template<typename BasisPtrType >
void createL2ProjectionStruct (const BasisPtrType cellBasis, const ordinal_type targetCubDegree)
 Initialize the ProjectionStruct for L2 projections. More...
 
template<typename BasisPtrType >
void createHGradProjectionStruct (const BasisPtrType cellBasis, const ordinal_type targetCubDegree, const ordinal_type targetGradCubDegre)
 Initialize the ProjectionStruct for HGRAD projections. More...
 
template<typename BasisPtrType >
void createHCurlProjectionStruct (const BasisPtrType cellBasis, const ordinal_type targetCubDegree, const ordinal_type targetCurlCubDegre)
 Initialize the ProjectionStruct for HCURL projections. More...
 
template<typename BasisPtrType >
void createHDivProjectionStruct (const BasisPtrType cellBasis, const ordinal_type targetCubDegree, const ordinal_type targetDivCubDegre)
 Initialize the ProjectionStruct for HDIV projections. More...
 
template<typename BasisPtrType >
void createHVolProjectionStruct (const BasisPtrType cellBasis, const ordinal_type targetCubDegree)
 Initialize the ProjectionStruct for HVOL (local-L2) projection. More...
 

Public Attributes

key_tag subCellTopologyKey
 
range_tag basisPointsRange
 
range_tag basisDerivPointsRange
 
range_tag targetPointsRange
 
range_tag targetDerivPointsRange
 
view_tag basisCubPoints
 
view_tag basisCubWeights
 
view_tag basisDerivCubPoints
 
view_tag basisDerivCubWeights
 
view_tag targetCubPoints
 
view_tag targetCubWeights
 
view_tag targetDerivCubPoints
 
view_tag targetDerivCubWeights
 
ordinal_type numBasisEvalPoints
 
ordinal_type numBasisDerivEvalPoints
 
ordinal_type numTargetEvalPoints
 
ordinal_type numTargetDerivEvalPoints
 
ordinal_type maxNumBasisEvalPoints
 
ordinal_type maxNumTargetEvalPoints
 
ordinal_type maxNumBasisDerivEvalPoints
 
ordinal_type maxNumTargetDerivEvalPoints
 

Static Public Attributes

static constexpr int numberSubCellDims = 4
 
static constexpr int maxSubCellsCount = 12
 

Detailed Description

template<typename SpT, typename ValueType>
class Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >

An helper class to compute the evaluation points and weights needed for performing projections.

In order to perform projections, the basis functions and the target function need to be evaluated at several sets of evaluation points (cubature points) defined on subcell entities (edges, faces, volumes). Depending on the projection, the evaluation of derivatives of the basis functions and of the target function may be needed as well. This class provides a struct to store the evaluation points/weights on the reference cell.

Use: create the proper ProjectionStruct rule by calling one of the functions: createL2ProjectionStruct, createHGradProjectionStruct, createHCurlProjectionStruct, createHDivProjectionStruct, createHVolProjectionStruct, depending on the type of projection wanted.

The created class is then used with the Projection Tools. See ProjectionTools class for more info.

Definition at line 87 of file Intrepid2_ProjectionStruct.hpp.

Member Function Documentation

template<typename SpT , typename ValueType >
template<typename BasisPtrType >
void Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::createHCurlProjectionStruct ( const BasisPtrType  cellBasis,
const ordinal_type  targetCubDegree,
const ordinal_type  targetCurlCubDegre 
)

Initialize the ProjectionStruct for HCURL projections.

Parameters
cellBasis[in] - HCURL basis functions for the projection
targetCubDegree[in] - degree of the cubature needed to integrate the target function
targetGradCubDegre[in] - degree of the cubature needed to integrate the derivative of target function

Definition at line 313 of file Intrepid2_ProjectionStructDef.hpp.

References Intrepid2::DefaultCubatureFactory::create().

template<typename SpT , typename ValueType >
template<typename BasisPtrType >
void Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::createHDivProjectionStruct ( const BasisPtrType  cellBasis,
const ordinal_type  targetCubDegree,
const ordinal_type  targetDivCubDegre 
)

Initialize the ProjectionStruct for HDIV projections.

Parameters
cellBasis[in] - HDIV basis functions for the projection
targetCubDegree[in] - degree of the cubature needed to integrate the target function
targetGradCubDegre[in] - degree of the cubature needed to integrate the derivative of target function

Definition at line 442 of file Intrepid2_ProjectionStructDef.hpp.

References Intrepid2::DefaultCubatureFactory::create(), and Intrepid2::Basis< ExecSpaceType, outputValueType, pointValueType >::getDofCount().

template<typename SpT , typename ValueType >
template<typename BasisPtrType >
void Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::createHGradProjectionStruct ( const BasisPtrType  cellBasis,
const ordinal_type  targetCubDegree,
const ordinal_type  targetGradCubDegre 
)

Initialize the ProjectionStruct for HGRAD projections.

Parameters
cellBasis[in] - HGRAD basis functions for the projection
targetCubDegree[in] - degree of the cubature needed to integrate the target function
targetGradCubDegre[in] - degree of the cubature needed to integrate the derivative of target function

Definition at line 196 of file Intrepid2_ProjectionStructDef.hpp.

References Intrepid2::DefaultCubatureFactory::create(), and Intrepid2::CellTools< ExecSpaceType >::getReferenceVertex().

template<typename SpT , typename ValueType >
template<typename BasisPtrType >
void Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::createHVolProjectionStruct ( const BasisPtrType  cellBasis,
const ordinal_type  targetCubDegree 
)

Initialize the ProjectionStruct for HVOL (local-L2) projection.

Parameters
cellBasis[in] - HVOL basis functions for the projection
targetCubDegree[in] - degree of the cubature needed to integrate the target function

Definition at line 558 of file Intrepid2_ProjectionStructDef.hpp.

References Intrepid2::DefaultCubatureFactory::create().

template<typename SpT , typename ValueType >
template<typename BasisPtrType >
void Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::createL2ProjectionStruct ( const BasisPtrType  cellBasis,
const ordinal_type  targetCubDegree 
)

Initialize the ProjectionStruct for L2 projections.

Parameters
cellBasis[in] - basis functions for the projection
targetCubDegree[in] - degree of the cubature needed to integrate the target function

Definition at line 67 of file Intrepid2_ProjectionStructDef.hpp.

References Intrepid2::DefaultCubatureFactory::create(), and Intrepid2::CellTools< ExecSpaceType >::getReferenceVertex().

template<typename SpT, typename ValueType>
view_type Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::getBasisDerivEvalPoints ( const ordinal_type  subCellDim,
const ordinal_type  subCellId 
)
inline

Returns the evaluation points for basis derivatives on a subcell.

P - num. evaluation points
D - spatial dimension
Parameters
subCellDim[in] - dimension of the subcell
subCellId[in] - ordinal of the subcell defined by cell topology
Returns
a rank-2 view (P,D) containing the basis derivatives evaluation points on the selected subcell

Definition at line 183 of file Intrepid2_ProjectionStruct.hpp.

template<typename SpT, typename ValueType>
view_type Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::getBasisDerivEvalWeights ( const ordinal_type  subCellDim,
const ordinal_type  subCellId 
)
inline

Returns the basis derivatives evaluation weights on a subcell.

P - num. evaluation points
Parameters
subCellDim[in] - dimension of the subcell
subCellId[in] - ordinal of the subcell defined by cell topology
Returns
a rank-1 view (P) containing the basis derivatives evaluation weights on the selected subcell

Definition at line 293 of file Intrepid2_ProjectionStruct.hpp.

template<typename SpT, typename ValueType>
const range_tag Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::getBasisDerivPointsRange ( ) const
inline

Returns the range tag of the derivative evaluation points on subcell.

Returns
the range tag of the basis derivative evaluation points corresponding on subcell

Definition at line 356 of file Intrepid2_ProjectionStruct.hpp.

template<typename SpT, typename ValueType>
view_type Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::getBasisEvalPoints ( const ordinal_type  subCellDim,
const ordinal_type  subCellId 
)
inline

Returns the basis evaluation points on a subcell.

P - num. evaluation points
D - spatial dimension
Parameters
subCellDim[in] - dimension of the subcell
subCellId[in] - ordinal of the subcell defined by cell topology
Returns
a rank-2 view (P,D) containing the basis evaluation points on the selected subcell

Definition at line 166 of file Intrepid2_ProjectionStruct.hpp.

template<typename SpT, typename ValueType>
view_type Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::getBasisEvalWeights ( const ordinal_type  subCellDim,
const ordinal_type  subCellId 
)
inline

Returns the basis evaluation weights on a subcell.

P - num. evaluation points
Parameters
subCellDim[in] - dimension of the subcell
subCellId[in] - ordinal of the subcell defined by cell topology
Returns
a rank-1 view (P) containing the basis evaluation weights on the selected subcell

Definition at line 277 of file Intrepid2_ProjectionStruct.hpp.

template<typename SpT, typename ValueType>
const range_tag Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::getBasisPointsRange ( ) const
inline

Returns the range tag of the basis evaluation points subcells.

Returns
the range tag of the basis evaluation points on subcells

Definition at line 334 of file Intrepid2_ProjectionStruct.hpp.

template<typename SpT, typename ValueType>
view_type Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::getDerivEvalPoints ( const ordinal_type  subCellDim,
const ordinal_type  subCellId,
EvalPointsType  type 
) const
inline

Returns the evaluation points for basis/target derivatives on a subcell.

P - num. evaluation points
D - spatial dimension
Parameters
subCellDim[in] - dimension of the subcell
subCellId[in] - ordinal of the subcell defined by cell topology
evalPointType[in] - enum selecting whether the points should be computed for the basis functions or for the target function
Returns
a rank-2 view (P,D) containing the basis/target derivatives evaluation points on the selected subcell

Definition at line 257 of file Intrepid2_ProjectionStruct.hpp.

template<typename SpT, typename ValueType>
const range_tag Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::getDerivPointsRange ( const EvalPointsType  type) const
inline

Returns the range tag of the basis/target derivative evaluation points on subcells.

Parameters
evalPointType[in] - enum selecting whether the points should be computed for the basis functions or for the target function
Returns
the range tag of the basis/target derivative evaluation points on subcells

Definition at line 366 of file Intrepid2_ProjectionStruct.hpp.

template<typename SpT, typename ValueType>
view_type Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::getEvalPoints ( const ordinal_type  subCellDim,
const ordinal_type  subCellId,
EvalPointsType  type 
) const
inline

Returns the basis/target evaluation points on a subcell.

P - num. evaluation points
D - spatial dimension
Parameters
subCellDim[in] - dimension of the subcell
subCellId[in] - ordinal of the subcell defined by cell topology
evalPointType[in] - enum selecting whether the points should be computed for the basis functions or for the target function
Returns
a rank-2 view (P,D) containing the basis/target function evaluation points on the selected subcell

Definition at line 236 of file Intrepid2_ProjectionStruct.hpp.

template<typename SpT, typename ValueType>
ordinal_type Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::getMaxNumDerivPoints ( const EvalPointsType  type) const
inline

Returns the maximum number of derivative evaluation points across all the subcells.

Parameters
evalPointType[in] - enum selecting whether the points should be computed for the basis functions or for the target function
Returns
the maximum number of the derivative evaluation points across all the subcells

Definition at line 134 of file Intrepid2_ProjectionStruct.hpp.

template<typename SpT, typename ValueType>
ordinal_type Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::getMaxNumEvalPoints ( const EvalPointsType  type) const
inline

Returns the maximum number of evaluation points across all the subcells.

Parameters
evalPointType[in] - enum selecting whether the points should be computed for the basis functions or for the target function
Returns
the maximum number of the evaluation points across all the subcells

Definition at line 147 of file Intrepid2_ProjectionStruct.hpp.

template<typename SpT, typename ValueType>
const range_tag Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::getPointsRange ( const EvalPointsType  type) const
inline

Returns the range tag of the basis/target evaluation points in subcells.

Parameters
evalPointType[in] - enum selecting whether the points should be computed for the basis functions or for the target function
Returns
the range tag of the basis/target evaluation points on subcells

Definition at line 344 of file Intrepid2_ProjectionStruct.hpp.

template<typename SpT, typename ValueType>
view_type Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::getTargetDerivEvalPoints ( const ordinal_type  subCellDim,
const ordinal_type  subCellId 
)
inline

Returns the points where to evaluate the derivatives of the target function on a subcell.

P - num. evaluation points
D - spatial dimension
Parameters
subCellDim[in] - dimension of the subcell
subCellId[in] - ordinal of the subcell defined by cell topology
Returns
a rank-2 view (P,D) containing the target derivatives evaluation points on the selected subcell

Definition at line 217 of file Intrepid2_ProjectionStruct.hpp.

template<typename SpT, typename ValueType>
view_type Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::getTargetDerivEvalWeights ( const ordinal_type  subCellDim,
const ordinal_type  subCellId 
)
inline

Returns the function derivatives evaluation weights on a subcell.

P - num. evaluation points
Parameters
subCellDim[in] - dimension of the subcell
subCellId[in] - ordinal of the subcell defined by cell topology
Returns
a rank-1 view (P) containing the target derivatives evaluation weights on the selected subcell

Definition at line 325 of file Intrepid2_ProjectionStruct.hpp.

template<typename SpT, typename ValueType>
const range_tag Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::getTargetDerivPointsRange ( ) const
inline

Returns the range tag of the target function derivative evaluation points on subcells.

Returns
the range of the target function derivative evaluation points corresponding on subcells

Definition at line 387 of file Intrepid2_ProjectionStruct.hpp.

template<typename SpT, typename ValueType>
view_type Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::getTargetEvalPoints ( const ordinal_type  subCellDim,
const ordinal_type  subCellId 
)
inline

Returns the points where to evaluate the target function on a subcell.

P - num. evaluation points
D - spatial dimension
Parameters
subCellDim[in] - dimension of the subcell
subCellId[in] - ordinal of the subcell defined by cell topology
Returns
a rank-2 view (P,D) containing the target evaluation points on the selected subcell

Definition at line 200 of file Intrepid2_ProjectionStruct.hpp.

template<typename SpT, typename ValueType>
view_type Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::getTargetEvalWeights ( const ordinal_type  subCellDim,
const ordinal_type  subCellId 
)
inline

Returns the function evaluation weights on a subcell.

P - num. evaluation points
Parameters
subCellDim[in] - dimension of the subcell
subCellId[in] - ordinal of the subcell defined by cell topology
Returns
a rank-1 view (P) containing the target evaluation weights on the selected subcell

Definition at line 309 of file Intrepid2_ProjectionStruct.hpp.

template<typename SpT, typename ValueType>
const range_tag Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::getTargetPointsRange ( ) const
inline

Returns the range of the target function evaluation points on subcells.

Returns
the range of the target function evaluation points corresponding on subcells

Definition at line 378 of file Intrepid2_ProjectionStruct.hpp.

template<typename SpT, typename ValueType>
const key_tag Intrepid2::Experimental::ProjectionStruct< SpT, ValueType >::getTopologyKey ( ) const
inline

Returns the key tag for subcells.

Returns
the key tag of the selected subcells

Definition at line 395 of file Intrepid2_ProjectionStruct.hpp.


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