Intrepid2
|
Tools to compute orientations for degrees-of-freedom. More...
#include <Intrepid2_OrientationTools.hpp>
Public Types | |
typedef Kokkos::View< double ****, DeviceType > | CoeffMatrixDataViewType |
subcell ordinal, orientation, matrix m x n | |
Public Member Functions | |
template<typename BasisHostType > | |
OrientationTools< DT > ::CoeffMatrixDataViewType | createCoeffMatrixInternal (const BasisHostType *basis, const bool inverse) |
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) |
template<typename BasisType > | |
OrientationTools< DT > ::CoeffMatrixDataViewType | createCoeffMatrix (const BasisType *basis) |
template<typename BasisType > | |
OrientationTools< DT > ::CoeffMatrixDataViewType | createInvCoeffMatrix (const BasisType *basis) |
Static Public Member Functions | |
template<typename BasisType > | |
static CoeffMatrixDataViewType | createCoeffMatrix (const BasisType *basis) |
Create coefficient matrix. More... | |
template<typename BasisType > | |
static CoeffMatrixDataViewType | createInvCoeffMatrix (const BasisType *basis) |
Create inverse of coefficient matrix. More... | |
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. More... | |
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. More... | |
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. More... | |
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(). More... | |
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. More... | |
Static Public Attributes | |
static std::map< std::pair < std::string, ordinal_type > , CoeffMatrixDataViewType > | ortCoeffData |
key :: basis name, order, value :: matrix data view type | |
static std::map< std::pair < std::string, ordinal_type > , CoeffMatrixDataViewType > | 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. | |
Tools to compute orientations for degrees-of-freedom.
Definition at line 378 of file Intrepid2_OrientationTools.hpp.
|
inlinestatic |
Create coefficient matrix.
basis | [in] - basis type |
|
inlinestatic |
Create inverse of coefficient matrix.
basis | [in] - basis type |
|
inlinestatic |
Compute orientations of cells in a workset.
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 32 of file Intrepid2_OrientationToolsDefModifyBasis.hpp.
Referenced by Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::initializeOrientations().
|
inlinestatic |
Modify basis due to orientation.
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 250 of file Intrepid2_OrientationToolsDefModifyBasis.hpp.
References Intrepid2::RealSpaceTools< DeviceType >::clone().
Referenced by Intrepid2::ProjectionTools< DeviceType >::projectField().
|
inlinestatic |
Modify basis due to orientation, applying the inverse of the operator applied in modifyBasisByOrientation().
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 350 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().
|
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.
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 335 of file Intrepid2_OrientationToolsDefModifyBasis.hpp.
|
inlinestatic |
Modify an assembled (C,F1,F2) matrix according to orientation of the cells.
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 436 of file Intrepid2_OrientationToolsDefModifyBasis.hpp.
References Intrepid2::RealSpaceTools< DeviceType >::clone().