Intrepid2
Public Member Functions | Static Public Member Functions | List of all members
Intrepid2::Impl::OrientationTools Class Reference

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

#include <Intrepid2_OrientationTools.hpp>

Public Member Functions

template<typename VT >
KOKKOS_INLINE_FUNCTION void getModifiedLinePoint (VT &ot, const VT pt, const ordinal_type ort)
 
template<typename VT >
KOKKOS_INLINE_FUNCTION void getModifiedTrianglePoint (VT &ot0, VT &ot1, const VT pt0, const VT pt1, const ordinal_type ort)
 
template<typename VT >
KOKKOS_INLINE_FUNCTION void getModifiedQuadrilateralPoint (VT &ot0, VT &ot1, const VT pt0, const VT pt1, const ordinal_type ort)
 
template<typename outPointViewType >
void getJacobianOfOrientationMap (outPointViewType jacobian, const shards::CellTopology cellTopo, const ordinal_type cellOrt)
 
template<typename outPointViewType >
KOKKOS_INLINE_FUNCTION void getJacobianOfOrientationMap (outPointViewType jacobian, const unsigned cellTopoKey, const ordinal_type cellOrt)
 

Static Public Member Functions

template<typename ValueType >
static KOKKOS_INLINE_FUNCTION void getModifiedLinePoint (ValueType &ot, const ValueType pt, const ordinal_type ort)
 Computes modified point for line segment. More...
 
template<typename ValueType >
static KOKKOS_INLINE_FUNCTION void getModifiedTrianglePoint (ValueType &ot0, ValueType &ot1, const ValueType pt0, const ValueType pt1, const ordinal_type ort)
 Computes modified point for triangle. More...
 
template<typename ValueType >
static KOKKOS_INLINE_FUNCTION void getModifiedQuadrilateralPoint (ValueType &ot0, ValueType &ot1, const ValueType pt0, const ValueType pt1, const ordinal_type ort)
 Computes modified point for quadrilateral. More...
 
template<typename outPointViewType , typename refPointViewType >
static void mapToModifiedReference (outPointViewType outPoints, const refPointViewType refPoints, const shards::CellTopology cellTopo, const ordinal_type cellOrt=0)
 Computes modified parameterization maps of 1- and 2-subcells with orientation. More...
 
template<typename outPointViewType , typename refPointViewType >
static KOKKOS_INLINE_FUNCTION void mapToModifiedReference (outPointViewType outPoints, const refPointViewType refPoints, const unsigned cellTopoKey, const ordinal_type cellOrt=0)
 Computes modified parameterization maps of 1- and 2-subcells with orientation. More...
 
template<typename JacobianViewType >
static KOKKOS_INLINE_FUNCTION void getLineJacobian (JacobianViewType jacobian, const ordinal_type ort)
 Computes Jacobian of orientation map for line segment. More...
 
template<typename JacobianViewType >
static KOKKOS_INLINE_FUNCTION void getTriangleJacobian (JacobianViewType jacobian, const ordinal_type ort)
 Computes Jacobian of orientation map for triangle. More...
 
template<typename JacobianViewType >
static KOKKOS_INLINE_FUNCTION void getQuadrilateralJacobian (JacobianViewType jacobian, const ordinal_type ort)
 Computes Jacobian of orientation map for quadrilateral. More...
 
template<typename JacobianViewType >
static void getJacobianOfOrientationMap (JacobianViewType jacobian, const shards::CellTopology cellTopo, const ordinal_type cellOrt)
 Computes jacobian of the parameterization maps of 1- and 2-subcells with orientation. More...
 
template<typename JacobianViewType >
static KOKKOS_INLINE_FUNCTION void getJacobianOfOrientationMap (JacobianViewType jacobian, const unsigned cellTopoKey, const ordinal_type cellOrt)
 Computes jacobian of the parameterization maps of 1- and 2-subcells with orientation. More...
 
template<typename TanViewType , typename ParamViewType >
static KOKKOS_INLINE_FUNCTION void getRefSubcellTangents (TanViewType tangents, const ParamViewType subCellParametrization, const unsigned subcellTopoKey, const ordinal_type subCellOrd, const ordinal_type ort)
 Computes the (oriented) subCell tangents. More...
 
template<typename TanNormViewType , typename ParamViewType >
static KOKKOS_INLINE_FUNCTION void getRefSideTangentsAndNormal (TanNormViewType tangentsAndNormal, const ParamViewType subCellParametrization, const unsigned subcellTopoKey, const ordinal_type subCellOrd, const ordinal_type ort)
 Computes the (oriented) tangents and normal of the side subCell The normals are only defined for sides (subCells of dimension D-1) and it requires the computation of tangents. More...
 
template<typename coordsViewType , typename subcellCoordsViewType , typename ParamViewType >
static KOKKOS_INLINE_FUNCTION void mapSubcellCoordsToRefCell (coordsViewType cellCoords, const subcellCoordsViewType subCellCoords, const ParamViewType subcellParametrization, const unsigned subcellTopoKey, const ordinal_type subCellOrd, const ordinal_type ort)
 Maps points defined on the subCell manifold into the parent Cell accounting for orientation. More...
 
template<typename OutputViewType , typename subcellBasisHostType , typename cellBasisHostType >
static void getCoeffMatrix_HGRAD (OutputViewType &output, const subcellBasisHostType &subcellBasis, const cellBasisHostType &cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt, const bool inverse=false)
 Compute orientation matrix for HGRAD basis for a given subcell and its reference basis. More...
 
template<typename OutputViewType , typename subcellBasisHostType , typename cellBasisHostType >
static void getCoeffMatrix_HCURL (OutputViewType &output, const subcellBasisHostType &subcellBasis, const cellBasisHostType &cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt, const bool inverse=false)
 Compute orientation matrix for HCURL basis for a given subcell and its reference basis. More...
 
template<typename OutputViewType , typename subcellBasisHostType , typename cellBasisHostType >
static void getCoeffMatrix_HDIV (OutputViewType &output, const subcellBasisHostType &subcellBasis, const cellBasisHostType &cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt, const bool inverse=false)
 Compute orientation matrix for HDIV basis for a given subcell and its reference basis. More...
 
template<typename OutputViewType , typename cellBasisHostType >
static void getCoeffMatrix_HVOL (OutputViewType &output, const cellBasisHostType &cellBasis, const ordinal_type cellOrt, const bool inverse=false)
 Compute orientation matrix for HVOL basis for a given (2D or 1D) cell and its reference basis. This method is used only for side basis. More...
 

Detailed Description

Tools to compute orientations for degrees-of-freedom.

Definition at line 114 of file Intrepid2_OrientationTools.hpp.

Member Function Documentation

template<typename OutputViewType , typename subcellBasisHostType , typename cellBasisHostType >
static void Intrepid2::Impl::OrientationTools::getCoeffMatrix_HCURL ( OutputViewType &  output,
const subcellBasisHostType &  subcellBasis,
const cellBasisHostType &  cellBasis,
const ordinal_type  subcellId,
const ordinal_type  subcellOrt,
const bool  inverse = false 
)
inlinestatic

Compute orientation matrix for HCURL basis for a given subcell and its reference basis.

Parameters
output[out] - rank 2 coefficient matrix
subcellBasis[in] - reference subcell basis function
cellBasis[in] - cell basis function
subcellId[in] - subcell Id in the cell topology
subcellOrt[in] - orientation number between 0 and 1
inverse[in] - boolean, when true the inverse of the orintation matrix is computed
template<typename OutputViewType , typename subcellBasisHostType , typename cellBasisHostType >
void Intrepid2::Impl::OrientationTools::getCoeffMatrix_HDIV ( OutputViewType &  output,
const subcellBasisHostType &  subcellBasis,
const cellBasisHostType &  cellBasis,
const ordinal_type  subcellId,
const ordinal_type  subcellOrt,
const bool  inverse = false 
)
inlinestatic

Compute orientation matrix for HDIV basis for a given subcell and its reference basis.

Parameters
output[out] - rank 2 orientation matrix
subcellBasis[in] - reference subcell basis
cellBasis[in] - cell basis function
subcellId[in] - subcell Id in the cell topology
subcellOrt[in] - orientation number between 0 and 1
inverse[in] - boolean, when true the inverse of the orintation matrix is computed

Definition at line 175 of file Intrepid2_OrientationToolsDefCoeffMatrix_HDIV.hpp.

References Intrepid2::RefSubcellParametrization< DeviceType >::get(), Intrepid2::PointTools::getLattice(), Intrepid2::PointTools::getLatticeSize(), getRefSideTangentsAndNormal(), and mapSubcellCoordsToRefCell().

template<typename OutputViewType , typename subcellBasisHostType , typename cellBasisHostType >
void Intrepid2::Impl::OrientationTools::getCoeffMatrix_HGRAD ( OutputViewType &  output,
const subcellBasisHostType &  subcellBasis,
const cellBasisHostType &  cellBasis,
const ordinal_type  subcellId,
const ordinal_type  subcellOrt,
const bool  inverse = false 
)
inlinestatic

Compute orientation matrix for HGRAD basis for a given subcell and its reference basis.

Parameters
output[out] - rank 2 coefficient matrix
subcellBasis[in] - reference subcell basis function
cellBasis[in] - cell basis function
subcellId[in] - subcell Id in the cell topology
subcellOrt[in] - orientation number between 0 and 1
inverse[in] - boolean, when true the inverse of the orintation matrix is computed
Parameters
subcellBasisthis is device view
cellBasisthis must be host basis object
subcellIdthis also must be host basis object

Definition at line 135 of file Intrepid2_OrientationToolsDefCoeffMatrix_HGRAD.hpp.

References Intrepid2::RefSubcellParametrization< DeviceType >::get(), Intrepid2::PointTools::getLattice(), Intrepid2::PointTools::getLatticeSize(), mapSubcellCoordsToRefCell(), and mapToModifiedReference().

template<typename OutputViewType , typename cellBasisHostType >
void Intrepid2::Impl::OrientationTools::getCoeffMatrix_HVOL ( OutputViewType &  output,
const cellBasisHostType &  cellBasis,
const ordinal_type  cellOrt,
const bool  inverse = false 
)
inlinestatic

Compute orientation matrix for HVOL basis for a given (2D or 1D) cell and its reference basis. This method is used only for side basis.

Parameters
output[out] - rank 2 coefficient matrix
cellBasis[in] - cell basis function
subcellOrt[in] - orientation number between 0 and 1
inverse[in] - boolean, when true the inverse of the orintation matrix is computed
Parameters
cellBasisthis is device view
cellOrtthis also must be host basis object

Definition at line 114 of file Intrepid2_OrientationToolsDefCoeffMatrix_HVOL.hpp.

References getJacobianOfOrientationMap(), Intrepid2::PointTools::getLattice(), and mapToModifiedReference().

template<typename JacobianViewType >
static void Intrepid2::Impl::OrientationTools::getJacobianOfOrientationMap ( JacobianViewType  jacobian,
const shards::CellTopology  cellTopo,
const ordinal_type  cellOrt 
)
inlinestatic

Computes jacobian of the parameterization maps of 1- and 2-subcells with orientation.

Parameters
jacobian[out] - rank-2 (D,D) array with jacobian of the orientation map
cellTopo[in] - cell topology of the parameterized domain (1- and 2-subcells)
cellOrt[in] - cell orientation number (zero is aligned with shards default configuration

Referenced by getCoeffMatrix_HVOL(), and getRefSubcellTangents().

template<typename JacobianViewType >
static KOKKOS_INLINE_FUNCTION void Intrepid2::Impl::OrientationTools::getJacobianOfOrientationMap ( JacobianViewType  jacobian,
const unsigned  cellTopoKey,
const ordinal_type  cellOrt 
)
static

Computes jacobian of the parameterization maps of 1- and 2-subcells with orientation.

Parameters
jacobian[out] - rank-2 (D,D) array with jacobian of the orientation map
cellTopoKey[in] - key of the cell topology of the parameterized domain (1- and 2-subcells)
cellOrt[in] - cell orientation number (zero is aligned with shards default configuration
template<typename JacobianViewType >
KOKKOS_INLINE_FUNCTION void Intrepid2::Impl::OrientationTools::getLineJacobian ( JacobianViewType  jacobian,
const ordinal_type  ort 
)
static

Computes Jacobian of orientation map for line segment.

Parameters
jacobian[out] - rank-2 (D,D) array with jacobian of the orientation map
ort[in] - orientation number between 0 and 1

Definition at line 91 of file Intrepid2_OrientationToolsDefModifyPoints.hpp.

template<typename ValueType >
static KOKKOS_INLINE_FUNCTION void Intrepid2::Impl::OrientationTools::getModifiedLinePoint ( ValueType &  ot,
const ValueType  pt,
const ordinal_type  ort 
)
static

Computes modified point for line segment.

Parameters
ot[out] - modified point value
pt[in] - input point in [-1.0 , 1.0]
ort[in] - orientation number between 0 and 1

Referenced by mapToModifiedReference().

template<typename ValueType >
static KOKKOS_INLINE_FUNCTION void Intrepid2::Impl::OrientationTools::getModifiedQuadrilateralPoint ( ValueType &  ot0,
ValueType &  ot1,
const ValueType  pt0,
const ValueType  pt1,
const ordinal_type  ort 
)
static

Computes modified point for quadrilateral.

Parameters
ot0[out] - modified coordinate 0
ot1[out] - modified coordinate 1
pt0[out] - input coordinate 0
pt1[out] - input coordinate 1
ort[in] - orientation number between 0 and 7

Referenced by mapToModifiedReference().

template<typename ValueType >
static KOKKOS_INLINE_FUNCTION void Intrepid2::Impl::OrientationTools::getModifiedTrianglePoint ( ValueType &  ot0,
ValueType &  ot1,
const ValueType  pt0,
const ValueType  pt1,
const ordinal_type  ort 
)
static

Computes modified point for triangle.

Parameters
ot0[out] - modified coordinate 0
ot1[out] - modified coordinate 1
pt0[out] - input coordinate 0
pt1[out] - input coordinate 1
ort[in] - orientation number between 0 and 5

Referenced by mapToModifiedReference().

template<typename JacobianViewType >
KOKKOS_INLINE_FUNCTION void Intrepid2::Impl::OrientationTools::getQuadrilateralJacobian ( JacobianViewType  jacobian,
const ordinal_type  ort 
)
static

Computes Jacobian of orientation map for quadrilateral.

Parameters
jacobian[out] - rank-2 (D,D) array with jacobian of the orientation map
ort[in] - orientation number between 0 and 7

Definition at line 220 of file Intrepid2_OrientationToolsDefModifyPoints.hpp.

template<typename TanNormViewType , typename ParamViewType >
KOKKOS_INLINE_FUNCTION void Intrepid2::Impl::OrientationTools::getRefSideTangentsAndNormal ( TanNormViewType  tangentsAndNormal,
const ParamViewType  subCellParametrization,
const unsigned  subcellTopoKey,
const ordinal_type  subCellOrd,
const ordinal_type  ort 
)
static

Computes the (oriented) tangents and normal of the side subCell The normals are only defined for sides (subCells of dimension D-1) and it requires the computation of tangents.

Parameters
tangentsAndNormal[out] - rank-2 (D,D), (D-1) tangents and 1 normal. D: parent cell dimension
subcellParam[in] - rank-2 (N, scD+1), subCells parametrization. N:number of subcells, scD: subCell dimension
cellTopoKey[in] - key of the cell topology of the parameterized domain (1- and 2-subcells)
subCellOrd[in] - position of the subCell among subCells of a given dimension
ort[in] - subCell orientation number

Definition at line 400 of file Intrepid2_OrientationToolsDefModifyPoints.hpp.

References getRefSubcellTangents().

Referenced by getCoeffMatrix_HDIV().

template<typename TanViewType , typename ParamViewType >
KOKKOS_INLINE_FUNCTION void Intrepid2::Impl::OrientationTools::getRefSubcellTangents ( TanViewType  tangents,
const ParamViewType  subCellParametrization,
const unsigned  subcellTopoKey,
const ordinal_type  subCellOrd,
const ordinal_type  ort 
)
static

Computes the (oriented) subCell tangents.

Parameters
tangents[out] - rank-2 (scD,D), tangents of the subcell. scD: subCell dimension, D: parent cell dimension
subcellParam[in] - rank-2 (N, scD+1), subCells parametrization. N:number of subcells, scD: subCell dimension
cellTopoKey[in] - key of the cell topology of the parameterized domain (1- and 2-subcells)
subCellOrd[in] - position of the subCell among subCells of a given dimension
ort[in] - subCell orientation number

Definition at line 376 of file Intrepid2_OrientationToolsDefModifyPoints.hpp.

References getJacobianOfOrientationMap().

Referenced by getRefSideTangentsAndNormal().

template<typename JacobianViewType >
KOKKOS_INLINE_FUNCTION void Intrepid2::Impl::OrientationTools::getTriangleJacobian ( JacobianViewType  jacobian,
const ordinal_type  ort 
)
static

Computes Jacobian of orientation map for triangle.

Parameters
jacobian[out] - rank-2 (D,D) array with jacobian of the orientation map
ort[in] - orientation number between 0 and 5

Definition at line 151 of file Intrepid2_OrientationToolsDefModifyPoints.hpp.

template<typename coordsViewType , typename subcellCoordsViewType , typename ParamViewType >
KOKKOS_INLINE_FUNCTION void Intrepid2::Impl::OrientationTools::mapSubcellCoordsToRefCell ( coordsViewType  cellCoords,
const subcellCoordsViewType  subCellCoords,
const ParamViewType  subcellParametrization,
const unsigned  subcellTopoKey,
const ordinal_type  subCellOrd,
const ordinal_type  ort 
)
static

Maps points defined on the subCell manifold into the parent Cell accounting for orientation.

Parameters
cellCoords[out] - rank-2 (P, D), P points mapped in the parent cell manifold.
subCellCoords[in] - rank-2 (P, scD), P points defined on the subCell manifold, scD: subCell dimension
subcellParam[in] - rank-2 (N, scD+1), subCells parametrization. N:number of subCells, scD: subCell dimension
cellTopoKey[in] - key of the cell topology of the parameterized domain (1- and 2-subCells)
subCellOrd[in] - position of the subCell among subCells of a given dimension
ort[in] - subCell orientation number

Definition at line 425 of file Intrepid2_OrientationToolsDefModifyPoints.hpp.

References mapToModifiedReference().

Referenced by getCoeffMatrix_HDIV(), and getCoeffMatrix_HGRAD().

template<typename outPointViewType , typename refPointViewType >
void Intrepid2::Impl::OrientationTools::mapToModifiedReference ( outPointViewType  outPoints,
const refPointViewType  refPoints,
const shards::CellTopology  cellTopo,
const ordinal_type  cellOrt = 0 
)
inlinestatic

Computes modified parameterization maps of 1- and 2-subcells with orientation.

Parameters
outPoints[out] - rank-2 (P,D2) array with points in 1D or 2D modified domain with orientation
refPoints[in] - rank-2 (P,D2) array with points in 1D or 2D parameter domain
cellTopo[in] - cell topology of the parameterized domain (1- and 2-subcells)
cellOrt[in] - cell orientation number (zero is aligned with shards default configuration

Definition at line 255 of file Intrepid2_OrientationToolsDefModifyPoints.hpp.

Referenced by getCoeffMatrix_HGRAD(), getCoeffMatrix_HVOL(), and mapSubcellCoordsToRefCell().

template<typename outPointViewType , typename refPointViewType >
KOKKOS_INLINE_FUNCTION void Intrepid2::Impl::OrientationTools::mapToModifiedReference ( outPointViewType  outPoints,
const refPointViewType  refPoints,
const unsigned  cellTopoKey,
const ordinal_type  cellOrt = 0 
)
static

Computes modified parameterization maps of 1- and 2-subcells with orientation.

Parameters
outPoints[out] - rank-2 (P,D2) array with points in 1D or 2D modified domain with orientation
refPoints[in] - rank-2 (P,D2) array with points in 1D or 2D parameter domain
cellTopoKey[in] - key of the cell topology of the parameterized domain (1- and 2-subcells)
cellOrt[in] - cell orientation number (zero is aligned with shards default configuration

Definition at line 281 of file Intrepid2_OrientationToolsDefModifyPoints.hpp.

References getModifiedLinePoint(), getModifiedQuadrilateralPoint(), and getModifiedTrianglePoint().


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