Intrepid2
Public Types | Public Member Functions | Static Public Member Functions | Static Public Attributes | Static Private Member Functions | List of all members
Intrepid2::OrientationTools< DeviceType > Class Template Reference

Tools to compute orientations for degrees-of-freedom. More...

#include <Intrepid2_OrientationTools.hpp>

Public Types

using CoeffMatrixDataViewType = Kokkos::View< double ****, DeviceType >
 subcell ordinal, orientation, matrix m x n
 
using KeyType = std::pair< const std::string, ordinal_type >
 key :: basis name, order, value :: matrix data view type
 
using OrtCoeffDataType = std::map< KeyType, Kokkos::View< double ****, DeviceType > >
 

Public Member Functions

template<typename BasisHostType >
void init_HGRAD (typename OrientationTools< DT >::CoeffMatrixDataViewType matData, BasisHostType const *cellBasis, const bool inverse)
 
template<typename BasisHostType >
void init_HCURL (typename OrientationTools< DT >::CoeffMatrixDataViewType matData, BasisHostType const *cellBasis, const bool inverse)
 
template<typename BasisHostType >
void init_HDIV (typename OrientationTools< DT >::CoeffMatrixDataViewType matData, BasisHostType const *cellBasis, const bool inverse)
 
template<typename BasisHostType >
void init_HVOL (typename OrientationTools< DT >::CoeffMatrixDataViewType matData, BasisHostType const *cellBasis, const bool inverse)
 

Static Public Member Functions

template<typename BasisType >
static CoeffMatrixDataViewType createCoeffMatrix (const BasisType *basis)
 Create coefficient matrix.
 
template<typename BasisType >
static CoeffMatrixDataViewType createInvCoeffMatrix (const BasisType *basis)
 Create inverse of coefficient matrix.
 
static void clearCoeffMatrix ()
 Clear coefficient matrix.
 
template<typename elemOrtValueType , class ... elemOrtProperties, typename elemNodeValueType , class ... elemNodeProperties>
static void getOrientation (Kokkos::DynRankView< elemOrtValueType, elemOrtProperties... > elemOrts, const Kokkos::DynRankView< elemNodeValueType, elemNodeProperties... > elemNodes, const shards::CellTopology cellTopo, bool isSide=false)
 Compute orientations of cells in a workset.
 
template<typename outputValueType , class ... outputProperties, typename inputValueType , class ... inputProperties, typename OrientationViewType , typename BasisType >
static void modifyBasisByOrientation (Kokkos::DynRankView< outputValueType, outputProperties... > output, const Kokkos::DynRankView< inputValueType, inputProperties... > input, const OrientationViewType orts, const BasisType *basis, const bool transpose=false)
 Modify basis due to orientation.
 
template<typename outputValueType , class ... outputProperties, typename inputValueType , class ... inputProperties, typename OrientationViewType , typename BasisType >
static void modifyBasisByOrientationTranspose (Kokkos::DynRankView< outputValueType, outputProperties... > output, const Kokkos::DynRankView< inputValueType, inputProperties... > input, const OrientationViewType orts, const BasisType *basis)
 Modify basis due to orientation, applying the transpose of the operator applied in modifyBasisByOrientation(). If the input provided represents basis coefficents in the global orientation, then this method will appropriately transform them to the local orientation.
 
template<typename outputValueType , class ... outputProperties, typename inputValueType , class ... inputProperties, typename OrientationViewType , typename BasisType >
static void modifyBasisByOrientationInverse (Kokkos::DynRankView< outputValueType, outputProperties... > output, const Kokkos::DynRankView< inputValueType, inputProperties... > input, const OrientationViewType orts, const BasisType *basis, const bool transpose=false)
 Modify basis due to orientation, applying the inverse of the operator applied in modifyBasisByOrientation().
 
template<typename outputValueType , class ... outputProperties, typename inputValueType , class ... inputProperties, typename OrientationViewType , typename BasisTypeLeft , typename BasisTypeRight >
static void modifyMatrixByOrientation (Kokkos::DynRankView< outputValueType, outputProperties... > output, const Kokkos::DynRankView< inputValueType, inputProperties... > input, const OrientationViewType orts, const BasisTypeLeft *basisLeft, const BasisTypeRight *basisRight)
 Modify an assembled (C,F1,F2) matrix according to orientation of the cells.
 

Static Public Attributes

static OrtCoeffDataType ortCoeffData
 
static OrtCoeffDataType ortInvCoeffData
 

Static Private Member Functions

template<typename BasisHostType >
static CoeffMatrixDataViewType createCoeffMatrixInternal (const BasisHostType *basis, const bool invTrans=false)
 
template<typename BasisHostType >
static void init_HGRAD (CoeffMatrixDataViewType matData, BasisHostType const *cellBasis, const bool inverse=false)
 Compute orientation matrix for HGRAD basis.
 
template<typename BasisHostType >
static void init_HCURL (CoeffMatrixDataViewType matData, BasisHostType const *cellBasis, const bool inverse=false)
 Compute orientation matrix for HCURL basis.
 
template<typename BasisHostType >
static void init_HDIV (CoeffMatrixDataViewType matData, BasisHostType const *cellBasis, const bool inverse=false)
 Compute orientation matrix for HDIV basis.
 
template<typename BasisHostType >
static void init_HVOL (CoeffMatrixDataViewType matData, BasisHostType const *cellBasis, const bool inverse=false)
 Compute orientation matrix for HVOL basis.
 

Detailed Description

template<typename DeviceType>
class Intrepid2::OrientationTools< DeviceType >

Tools to compute orientations for degrees-of-freedom.

Definition at line 378 of file Intrepid2_OrientationTools.hpp.

Member Typedef Documentation

◆ CoeffMatrixDataViewType

template<typename DeviceType >
using Intrepid2::OrientationTools< DeviceType >::CoeffMatrixDataViewType = Kokkos::View<double****,DeviceType>

subcell ordinal, orientation, matrix m x n

Definition at line 383 of file Intrepid2_OrientationTools.hpp.

◆ KeyType

template<typename DeviceType >
using Intrepid2::OrientationTools< DeviceType >::KeyType = std::pair<const std::string,ordinal_type>

key :: basis name, order, value :: matrix data view type

Definition at line 388 of file Intrepid2_OrientationTools.hpp.

◆ OrtCoeffDataType

template<typename DeviceType >
using Intrepid2::OrientationTools< DeviceType >::OrtCoeffDataType = std::map<KeyType,Kokkos::View<double****,DeviceType> >

Definition at line 389 of file Intrepid2_OrientationTools.hpp.

Member Function Documentation

◆ clearCoeffMatrix()

template<typename DT >
void Intrepid2::OrientationTools< DT >::clearCoeffMatrix ( )
inlinestatic

Clear coefficient matrix.

Definition at line 318 of file Intrepid2_OrientationToolsDefMatrixData.hpp.

◆ createCoeffMatrix()

template<typename DT >
template<typename BasisType >
OrientationTools< DT >::CoeffMatrixDataViewType Intrepid2::OrientationTools< DT >::createCoeffMatrix ( const BasisType *  basis)
inlinestatic

Create coefficient matrix.

Parameters
basis[in] - basis type

Definition at line 270 of file Intrepid2_OrientationToolsDefMatrixData.hpp.

◆ createCoeffMatrixInternal()

template<typename DT >
template<typename BasisHostType >
OrientationTools< DT >::CoeffMatrixDataViewType Intrepid2::OrientationTools< DT >::createCoeffMatrixInternal ( const BasisHostType *  basis,
const bool  invTrans = false 
)
inlinestaticprivate

Definition at line 28 of file Intrepid2_OrientationToolsDefMatrixData.hpp.

◆ createInvCoeffMatrix()

template<typename DT >
template<typename BasisType >
OrientationTools< DT >::CoeffMatrixDataViewType Intrepid2::OrientationTools< DT >::createInvCoeffMatrix ( const BasisType *  basis)
inlinestatic

Create inverse of coefficient matrix.

Parameters
basis[in] - basis type

Definition at line 295 of file Intrepid2_OrientationToolsDefMatrixData.hpp.

◆ getOrientation()

template<typename DT >
template<typename elemOrtValueType , class ... elemOrtProperties, typename elemNodeValueType , class ... elemNodeProperties>
void Intrepid2::OrientationTools< DT >::getOrientation ( Kokkos::DynRankView< elemOrtValueType, elemOrtProperties... >  elemOrts,
const Kokkos::DynRankView< elemNodeValueType, elemNodeProperties... >  elemNodes,
const shards::CellTopology  cellTopo,
bool  isSide = false 
)
inlinestatic

Compute orientations of cells in a workset.

Parameters
elemOrts[out] - cell orientations
elemNodes[in] - node coordinates
cellTopo[in] - shards cell topology
isSide[in] - boolean, whether the cell is a side

Definition at line 31 of file Intrepid2_OrientationToolsDefModifyBasis.hpp.

Referenced by Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::initializeOrientations().

◆ init_HCURL()

template<typename DeviceType >
template<typename BasisHostType >
void Intrepid2::OrientationTools< DeviceType >::init_HCURL ( typename OrientationTools< DT >::CoeffMatrixDataViewType  matData,
BasisHostType const *  cellBasis,
const bool  inverse 
)

Definition at line 136 of file Intrepid2_OrientationToolsDefMatrixData.hpp.

◆ init_HDIV()

template<typename DeviceType >
template<typename BasisHostType >
void Intrepid2::OrientationTools< DeviceType >::init_HDIV ( typename OrientationTools< DT >::CoeffMatrixDataViewType  matData,
BasisHostType const *  cellBasis,
const bool  inverse 
)

Definition at line 196 of file Intrepid2_OrientationToolsDefMatrixData.hpp.

◆ init_HGRAD()

template<typename DeviceType >
template<typename BasisHostType >
void Intrepid2::OrientationTools< DeviceType >::init_HGRAD ( typename OrientationTools< DT >::CoeffMatrixDataViewType  matData,
BasisHostType const *  cellBasis,
const bool  inverse 
)

Definition at line 72 of file Intrepid2_OrientationToolsDefMatrixData.hpp.

◆ init_HVOL()

template<typename DeviceType >
template<typename BasisHostType >
void Intrepid2::OrientationTools< DeviceType >::init_HVOL ( typename OrientationTools< DT >::CoeffMatrixDataViewType  matData,
BasisHostType const *  cellBasis,
const bool  inverse 
)

Definition at line 230 of file Intrepid2_OrientationToolsDefMatrixData.hpp.

◆ modifyBasisByOrientation()

template<typename DT >
template<typename outputValueType , class ... outputProperties, typename inputValueType , class ... inputProperties, typename OrientationViewType , typename BasisType >
void Intrepid2::OrientationTools< DT >::modifyBasisByOrientation ( Kokkos::DynRankView< outputValueType, outputProperties... >  output,
const Kokkos::DynRankView< inputValueType, inputProperties... >  input,
const OrientationViewType  orts,
const BasisType *  basis,
const bool  transpose = false 
)
inlinestatic

Modify basis due to orientation.

Parameters
output[out] - output array, of shape (C,F,P[,D])
input[in] - input array, of shape (C,F,P[,D]) or (F,P[,D])
orts[in] - orientations, of shape (C)
basis[in] - basis of cardinality F
transpose[in] - boolean, when true the transpose of the orintation matrix is applied

Definition at line 249 of file Intrepid2_OrientationToolsDefModifyBasis.hpp.

References Intrepid2::RealSpaceTools< DeviceType >::clone().

Referenced by Intrepid2::ProjectionTools< DeviceType >::projectField().

◆ modifyBasisByOrientationInverse()

template<typename DT >
template<typename outputValueType , class ... outputProperties, typename inputValueType , class ... inputProperties, typename OrientationViewType , typename BasisType >
void Intrepid2::OrientationTools< DT >::modifyBasisByOrientationInverse ( Kokkos::DynRankView< outputValueType, outputProperties... >  output,
const Kokkos::DynRankView< inputValueType, inputProperties... >  input,
const OrientationViewType  orts,
const BasisType *  basis,
const bool  transpose = false 
)
inlinestatic

Modify basis due to orientation, applying the inverse of the operator applied in modifyBasisByOrientation().

Parameters
output[out] - output array, of shape (C,F,P[,D])
input[in] - input array, of shape (C,F,P[,D]) or (F,P[,D])
orts[in] - orientations, of shape (C)
basis[in] - basis of cardinality F
transpose[in] - boolean, when true the transpose of the inverse of the operatore applied

Definition at line 349 of file Intrepid2_OrientationToolsDefModifyBasis.hpp.

References Intrepid2::RealSpaceTools< DeviceType >::clone().

Referenced by Intrepid2::LagrangianInterpolation< DeviceType >::getBasisCoeffs(), Intrepid2::ProjectionTools< DeviceType >::getHCurlBasisCoeffs(), Intrepid2::ProjectionTools< DeviceType >::getHDivBasisCoeffs(), Intrepid2::ProjectionTools< DeviceType >::getHGradBasisCoeffs(), Intrepid2::ProjectionTools< DeviceType >::getL2BasisCoeffs(), Intrepid2::ProjectionTools< DeviceType >::getL2DGBasisCoeffs(), and Intrepid2::LagrangianTools< DeviceType >::getOrientedDofCoeffs().

◆ modifyBasisByOrientationTranspose()

template<typename DT >
template<typename outputValueType , class ... outputProperties, typename inputValueType , class ... inputProperties, typename OrientationViewType , typename BasisType >
void Intrepid2::OrientationTools< DT >::modifyBasisByOrientationTranspose ( Kokkos::DynRankView< outputValueType, outputProperties... >  output,
const Kokkos::DynRankView< inputValueType, inputProperties... >  input,
const OrientationViewType  orts,
const BasisType *  basis 
)
inlinestatic

Modify basis due to orientation, applying the transpose of the operator applied in modifyBasisByOrientation(). If the input provided represents basis coefficents in the global orientation, then this method will appropriately transform them to the local orientation.

Parameters
output[out] - output array, of shape (C,F,P[,D])
input[in] - input array, of shape (C,F,P[,D]) or (F,P[,D])
orts[in] - orientations, of shape (C)
basis[in] - basis of cardinality F

Definition at line 334 of file Intrepid2_OrientationToolsDefModifyBasis.hpp.

◆ modifyMatrixByOrientation()

template<typename DT >
template<typename outputValueType , class ... outputProperties, typename inputValueType , class ... inputProperties, typename OrientationViewType , typename BasisTypeLeft , typename BasisTypeRight >
void Intrepid2::OrientationTools< DT >::modifyMatrixByOrientation ( Kokkos::DynRankView< outputValueType, outputProperties... >  output,
const Kokkos::DynRankView< inputValueType, inputProperties... >  input,
const OrientationViewType  orts,
const BasisTypeLeft *  basisLeft,
const BasisTypeRight *  basisRight 
)
inlinestatic

Modify an assembled (C,F1,F2) matrix according to orientation of the cells.

Parameters
output[out] - output array, shape (C,F1,F2)
input[in] - input array, shape (C,F1,F2)
orts[in] - orientations, shape (C)
basisLeft[in] - basis with cardinality F1
basisRight[in] - basis with cardinality F2

Definition at line 435 of file Intrepid2_OrientationToolsDefModifyBasis.hpp.

References Intrepid2::RealSpaceTools< DeviceType >::clone().

Member Data Documentation

◆ ortCoeffData

template<typename T >
OrientationTools< T >::OrtCoeffDataType Intrepid2::OrientationTools< T >::ortCoeffData
static

Definition at line 391 of file Intrepid2_OrientationTools.hpp.

◆ ortInvCoeffData

template<typename T >
OrientationTools< T >::OrtCoeffDataType Intrepid2::OrientationTools< T >::ortInvCoeffData
static

Definition at line 392 of file Intrepid2_OrientationTools.hpp.


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