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 OutputViewType , typename subcellBasisType , typename cellBasisType >
static void getCoeffMatrix_HGRAD (OutputViewType &output, const subcellBasisType subcellBasis, const cellBasisType cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt)
 Compute coefficient matrix for HGRAD by collocating point values. More...
 
template<typename OutputViewType , typename subcellBasisType , typename cellBasisType >
static void getCoeffMatrix_HCURL (OutputViewType &output, const subcellBasisType subcellBasis, const cellBasisType cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt)
 Compute coefficient matrix for HCURL by collocating point values. More...
 
template<typename OutputViewType , typename subcellBasisType , typename cellBasisType >
static void getCoeffMatrix_HDIV (OutputViewType &output, const subcellBasisType subcellBasis, const cellBasisType cellBasis, const ordinal_type subcellId, const ordinal_type subcellOrt)
 Compute coefficient matrix for HDIV by collocating point values. 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 subcellBasisType , typename cellBasisType >
void Intrepid2::Impl::OrientationTools::getCoeffMatrix_HCURL ( OutputViewType &  output,
const subcellBasisType  subcellBasis,
const cellBasisType  cellBasis,
const ordinal_type  subcellId,
const ordinal_type  subcellOrt 
)
inlinestatic

Compute coefficient matrix for HCURL by collocating point values.

Parameters
output[out] - rank 2 coefficient matrix
subcellBasis[in] - subcell basis function
cellBasis[in] - cell basis function
subcellId[in] - subcell Id in the cell topology
subcellOrt[in] - orientation number between 0 and 1

Definition at line 77 of file Intrepid2_OrientationToolsDefCoeffMatrix_HCURL.hpp.

References getJacobianOfOrientationMap(), Intrepid2::CellTools< ExecSpaceType >::getReferenceEdgeTangent(), Intrepid2::CellTools< ExecSpaceType >::getReferenceFaceTangents(), mapToModifiedReference(), and Intrepid2::CellTools< ExecSpaceType >::mapToReferenceSubcell().

template<typename OutputViewType , typename subcellBasisType , typename cellBasisType >
void Intrepid2::Impl::OrientationTools::getCoeffMatrix_HDIV ( OutputViewType &  output,
const subcellBasisType  subcellBasis,
const cellBasisType  cellBasis,
const ordinal_type  subcellId,
const ordinal_type  subcellOrt 
)
inlinestatic

Compute coefficient matrix for HDIV by collocating point values.

Parameters
output[out] - rank 2 coefficient matrix
subcellBasis[in] - subcell basis function
cellBasis[in] - cell basis function
subcellId[in] - subcell Id in the cell topology
subcellOrt[in] - orientation number between 0 and 1

Definition at line 81 of file Intrepid2_OrientationToolsDefCoeffMatrix_HDIV.hpp.

References getJacobianOfOrientationMap(), Intrepid2::CellTools< ExecSpaceType >::getReferenceSideNormal(), mapToModifiedReference(), and Intrepid2::CellTools< ExecSpaceType >::mapToReferenceSubcell().

template<typename OutputViewType , typename subcellBasisType , typename cellBasisType >
void Intrepid2::Impl::OrientationTools::getCoeffMatrix_HGRAD ( OutputViewType &  output,
const subcellBasisType  subcellBasis,
const cellBasisType  cellBasis,
const ordinal_type  subcellId,
const ordinal_type  subcellOrt 
)
inlinestatic

Compute coefficient matrix for HGRAD by collocating point values.

Parameters
output[out] - rank 2 coefficient matrix
subcellBasis[in] - subcell basis function
cellBasis[in] - cell basis function
subcellId[in] - subcell Id in the cell topology
subcellOrt[in] - orientation number between 0 and 1

Definition at line 73 of file Intrepid2_OrientationToolsDefCoeffMatrix_HGRAD.hpp.

References mapToModifiedReference(), and Intrepid2::CellTools< ExecSpaceType >::mapToReferenceSubcell().

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_HCURL(), and getCoeffMatrix_HDIV().

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 219 of file Intrepid2_OrientationToolsDefModifyPoints.hpp.

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 150 of file Intrepid2_OrientationToolsDefModifyPoints.hpp.

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 254 of file Intrepid2_OrientationToolsDefModifyPoints.hpp.

Referenced by getCoeffMatrix_HCURL(), getCoeffMatrix_HDIV(), and getCoeffMatrix_HGRAD().

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 280 of file Intrepid2_OrientationToolsDefModifyPoints.hpp.

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


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