Intrepid2
Static Public Member Functions | List of all members
Intrepid2::Impl::CellTools::Serial Struct Reference

Static Public Member Functions

template<typename jacobianViewType , typename basisGradViewType , typename nodeViewType >
static KOKKOS_INLINE_FUNCTION void computeJacobian (const jacobianViewType &jacobian, const basisGradViewType &grads, const nodeViewType &nodes)
 
template<typename PointViewType , typename basisValViewType , typename nodeViewType >
static KOKKOS_INLINE_FUNCTION void mapToPhysicalFrame (const PointViewType &point, const basisValViewType &vals, const nodeViewType &nodes)
 
template<typename implBasisType , typename refPointViewType , typename phyPointViewType , typename nodeViewType >
static KOKKOS_INLINE_FUNCTION bool mapToReferenceFrame (const refPointViewType &refPoint, const phyPointViewType &physPoint, const nodeViewType &nodes, const double tol=tolerence())
 Computation of $ F^{-1}_{c} $, the inverse of the reference-to-physical frame map. More...
 
template<typename refEdgeTanViewType , typename ParamViewType >
static KOKKOS_INLINE_FUNCTION void getReferenceEdgeTangent (const refEdgeTanViewType refEdgeTangent, const ParamViewType edgeParametrization, const ordinal_type edgeOrdinal)
 
template<typename refFaceTanViewType , typename ParamViewType >
static KOKKOS_INLINE_FUNCTION void getReferenceFaceTangent (const refFaceTanViewType refFaceTanU, const refFaceTanViewType refFaceTanV, const ParamViewType faceParametrization, const ordinal_type faceOrdinal)
 
template<typename edgeTangentViewType , typename ParamViewType , typename jacobianViewType >
static KOKKOS_INLINE_FUNCTION void getPhysicalEdgeTangent (const edgeTangentViewType edgeTangent, const ParamViewType edgeParametrization, const jacobianViewType jacobian, const ordinal_type edgeOrdinal)
 
template<typename faceTanViewType , typename ParamViewType , typename jacobianViewType >
static KOKKOS_INLINE_FUNCTION void getPhysicalFaceTangents (faceTanViewType faceTanU, faceTanViewType faceTanV, const ParamViewType faceParametrization, const jacobianViewType jacobian, const ordinal_type faceOrdinal)
 
template<typename faceNormalViewType , typename ParamViewType , typename jacobianViewType >
static KOKKOS_INLINE_FUNCTION void getPhysicalFaceNormal (const faceNormalViewType faceNormal, const ParamViewType faceParametrization, const jacobianViewType jacobian, const ordinal_type faceOrdinal)
 
template<typename sideNormalViewType , typename ParamViewType , typename jacobianViewType >
static KOKKOS_INLINE_FUNCTION void getPhysicalSideNormal (const sideNormalViewType sideNormal, const ParamViewType sideParametrization, const jacobianViewType jacobian, const ordinal_type sideOrdinal)
 

Detailed Description

Definition at line 72 of file Intrepid2_CellTools_Serial.hpp.

Member Function Documentation

template<typename implBasisType , typename refPointViewType , typename phyPointViewType , typename nodeViewType >
static KOKKOS_INLINE_FUNCTION bool Intrepid2::Impl::CellTools::Serial::mapToReferenceFrame ( const refPointViewType &  refPoint,
const phyPointViewType &  physPoint,
const nodeViewType &  nodes,
const double  tol = tolerence() 
)
inlinestatic

Computation of $ F^{-1}_{c} $, the inverse of the reference-to-physical frame map.

Applies $ F^{-1}_{c} $ for a cell workset to a single point (P,D) view. Computes a rank-2 (P,D) view such that

\[ \mbox{refPoints}(p,d) = \Big(F^{-1}_c(physPoint(p,*)) \Big)_d \]

Requires pointer to HGrad basis that defines reference to physical cell mapping. See Section Reference-to-physical cell mapping for definition of the mapping function. Note that the physical point can be in a space with larger dimension, e.g. if mapping a reference quadrilateral into a side in 3D space.

Warning
The array physPoints represents an arbitrary set (or sets) of points in the physical frame that are not required to belong in the physical cell (cells) that define(s) the reference to physical mapping. As a result, the images of these points in the reference frame are not necessarily contained in the reference cell corresponding to the specified cell topology.
Parameters
refPoint[out] - rank-2 view with dimension (refD) with the reference point
physPoint[in] - rank-2 view with dimension (physD) with the point in physical frame
cellNodes[in] - rank-2 array with dimensions (N,physD) with the nodes of the cell containing the physical point
Returns
a boolean set to true if the algorithm converged, to false otherise.

Definition at line 152 of file Intrepid2_CellTools_Serial.hpp.

References Intrepid2::Parameters::MaxNewton.


The documentation for this struct was generated from the following file: