Intrepid2
|
A stateless class for operations on cell data. Provides methods for: More...
#include <Intrepid2_CellTools.hpp>
Public Member Functions | |
CellTools ()=default | |
Default constructor. | |
~CellTools ()=default | |
Destructor. | |
template<typename InCellViewType , typename InputViewType> | |
void | checkPointwiseInclusion (InCellViewType inCell, const InputViewType points, const shards::CellTopology cellTopo, const typename ScalarTraits< typename InputViewType::value_type >::scalar_type threshold) |
template<typename Scalar > | |
void | setJacobianDividedByDet (Data< Scalar, DeviceType > &jacobianDividedByDet, const Data< Scalar, DeviceType > &jacobian, const Data< Scalar, DeviceType > &jacobianDetInv) |
template<typename PhysPointViewType , typename RefPointViewType , typename WorksetType , typename HGradBasisPtrType > | |
void | mapToPhysicalFrame (PhysPointViewType physPoints, const RefPointViewType refPoints, const WorksetType worksetCell, const HGradBasisPtrType basis) |
template<typename refSubcellViewType , typename paramViewType > | |
void | mapToReferenceSubcell (refSubcellViewType refSubcellPoints, const paramViewType paramPoints, const ordinal_type subcellDim, const ordinal_type subcellOrd, const shards::CellTopology parentCell) |
template<typename refSubcellViewType , typename paramViewType > | |
void | mapToReferenceSubcell (refSubcellViewType refSubcellPoints, const paramViewType paramPoints, const typename RefSubcellParametrization< DeviceType >::ConstViewType subcellParametrization, const ordinal_type subcellOrd) |
template<typename refSubcellViewType , typename paramViewType , typename ordViewType > | |
void | mapToReferenceSubcell (refSubcellViewType refSubcellPoints, const paramViewType paramPoints, const typename RefSubcellParametrization< DeviceType >::ConstViewType subcellParametrization, const ordViewType subcellOrd) |
Static Public Member Functions | |
static bool | hasReferenceCell (const shards::CellTopology cellTopo) |
Checks if a cell topology has reference cell. More... | |
template<typename JacobianViewType , typename PointViewType , typename WorksetType , typename HGradBasisType > | |
static void | setJacobian (JacobianViewType jacobian, const PointViewType points, const WorksetType worksetCell, const Teuchos::RCP< HGradBasisType > basis, const int startCell=0, const int endCell=-1) |
Computes the Jacobian matrix DF of the reference-to-physical frame map F. More... | |
template<typename JacobianViewType , typename BasisGradientsType , typename WorksetType > | |
static void | setJacobian (JacobianViewType jacobian, const WorksetType worksetCell, const BasisGradientsType gradients, const int startCell=0, const int endCell=-1) |
Computes the Jacobian matrix DF of the reference-to-physical frame map F. More... | |
template<typename JacobianViewType , typename PointViewType , typename WorksetCellViewType > | |
static void | setJacobian (JacobianViewType jacobian, const PointViewType points, const WorksetCellViewType worksetCell, const shards::CellTopology cellTopo) |
Computes the Jacobian matrix DF of the reference-to-physical frame map F. More... | |
template<typename JacobianInvViewType , typename JacobianViewType > | |
static void | setJacobianInv (JacobianInvViewType jacobianInv, const JacobianViewType jacobian) |
Computes the inverse of the Jacobian matrix DF of the reference-to-physical frame map F. More... | |
template<typename JacobianDetViewType , typename JacobianViewType > | |
static void | setJacobianDet (JacobianDetViewType jacobianDet, const JacobianViewType jacobian) |
Computes the determinant of the Jacobian matrix DF of the reference-to-physical frame map F. More... | |
template<class PointScalar > | |
static Data< PointScalar, DeviceType > | allocateJacobianDet (const Data< PointScalar, DeviceType > &jacobian) |
Allocates and returns a Data container suitable for storing determinants corresponding to the Jacobians in the Data container provided. More... | |
template<class PointScalar > | |
static Data< PointScalar, DeviceType > | allocateJacobianInv (const Data< PointScalar, DeviceType > &jacobian) |
Allocates and returns a Data container suitable for storing inverses corresponding to the Jacobians in the Data container provided. More... | |
template<class PointScalar > | |
static void | setJacobianDet (Data< PointScalar, DeviceType > &jacobianDet, const Data< PointScalar, DeviceType > &jacobian) |
Computes determinants corresponding to the Jacobians in the Data container provided. More... | |
template<class PointScalar > | |
static void | setJacobianDetInv (Data< PointScalar, DeviceType > &jacobianDetInv, const Data< PointScalar, DeviceType > &jacobian) |
Computes reciprocals of determinants corresponding to the Jacobians in the Data container provided. More... | |
template<class PointScalar > | |
static void | setJacobianInv (Data< PointScalar, DeviceType > &jacobianInv, const Data< PointScalar, DeviceType > &jacobian) |
Computes determinants corresponding to the Jacobians in the Data container provided. More... | |
template<class PointScalar > | |
static void | setJacobianDividedByDet (Data< PointScalar, DeviceType > &jacobianDividedByDet, const Data< PointScalar, DeviceType > &jacobian, const Data< PointScalar, DeviceType > &jacobianDetInv) |
Multiplies the Jacobian with shape (C,P,D,D) by the reciprocals of the determinants, with shape (C,P), entrywise. More... | |
template<typename cellCenterValueType , class... cellCenterProperties> | |
static void | getReferenceCellCenter (Kokkos::DynRankView< cellCenterValueType, cellCenterProperties...> cellCenter, const shards::CellTopology cell) |
Computes the Cartesian coordinates of reference cell barycenter. More... | |
template<typename cellVertexValueType , class... cellVertexProperties> | |
static void | getReferenceVertex (Kokkos::DynRankView< cellVertexValueType, cellVertexProperties...> cellVertex, const shards::CellTopology cell, const ordinal_type vertexOrd) |
Retrieves the Cartesian coordinates of a reference cell vertex. More... | |
template<typename subcellVertexValueType , class... subcellVertexProperties> | |
static void | getReferenceSubcellVertices (Kokkos::DynRankView< subcellVertexValueType, subcellVertexProperties...> subcellVertices, const ordinal_type subcellDim, const ordinal_type subcellOrd, const shards::CellTopology parentCell) |
Retrieves the Cartesian coordinates of all vertices of a reference subcell. More... | |
template<typename cellNodeValueType , class... cellNodeProperties> | |
static void | getReferenceNode (Kokkos::DynRankView< cellNodeValueType, cellNodeProperties...> cellNode, const shards::CellTopology cell, const ordinal_type nodeOrd) |
Retrieves the Cartesian coordinates of a reference cell node. More... | |
template<typename SubcellNodeViewType > | |
static void | getReferenceSubcellNodes (SubcellNodeViewType subcellNodes, const ordinal_type subcellDim, const ordinal_type subcellOrd, const shards::CellTopology parentCell) |
Retrieves the Cartesian coordinates of all nodes of a reference subcell. More... | |
template<typename RefEdgeTangentViewType > | |
static void | getReferenceEdgeTangent (RefEdgeTangentViewType refEdgeTangent, const ordinal_type edgeOrd, const shards::CellTopology parentCell) |
Computes constant tangent vectors to edges of 2D or 3D reference cells. More... | |
template<typename RefFaceTanViewType > | |
static void | getReferenceFaceTangents (RefFaceTanViewType refFaceTanU, RefFaceTanViewType refFaceTanV, const ordinal_type faceOrd, const shards::CellTopology parentCell) |
Computes pairs of constant tangent vectors to faces of a 3D reference cells. More... | |
template<typename RefSideNormalViewType > | |
static void | getReferenceSideNormal (RefSideNormalViewType refSideNormal, const ordinal_type sideOrd, const shards::CellTopology parentCell) |
Computes constant normal vectors to sides of 2D or 3D reference cells. More... | |
template<typename RefFaceNormalViewType > | |
static void | getReferenceFaceNormal (RefFaceNormalViewType refFaceNormal, const ordinal_type faceOrd, const shards::CellTopology parentCell) |
Computes constant normal vectors to faces of 3D reference cell. More... | |
template<typename edgeTangentValueType , class... edgeTangentProperties, typename worksetJacobianValueType , class... worksetJacobianProperties> | |
static void | getPhysicalEdgeTangents (Kokkos::DynRankView< edgeTangentValueType, edgeTangentProperties...> edgeTangents, const Kokkos::DynRankView< worksetJacobianValueType, worksetJacobianProperties...> worksetJacobians, const ordinal_type worksetEdgeOrd, const shards::CellTopology parentCell) |
Computes non-normalized tangent vectors to physical edges in an edge workset ; (see Subcell worksets for definition of edge worksets). More... | |
template<typename edgeTangentValueType , class... edgeTangentProperties, typename worksetJacobianValueType , class... worksetJacobianProperties, typename edgeOrdValueType , class... edgeOrdProperties> | |
static void | getPhysicalEdgeTangents (Kokkos::DynRankView< edgeTangentValueType, edgeTangentProperties...> edgeTangents, const Kokkos::DynRankView< worksetJacobianValueType, worksetJacobianProperties...> worksetJacobians, const Kokkos::DynRankView< edgeOrdValueType, edgeOrdProperties...> worksetEdgeOrds, const shards::CellTopology parentCell) |
Computes non-normalized tangent vectors to physical edges in an edge workset ; (see Subcell worksets for definition of edge worksets). More... | |
template<typename faceTanValueType , class... faceTanProperties, typename worksetJacobianValueType , class... worksetJacobianProperties> | |
static void | getPhysicalFaceTangents (Kokkos::DynRankView< faceTanValueType, faceTanProperties...> faceTanU, Kokkos::DynRankView< faceTanValueType, faceTanProperties...> faceTanV, const Kokkos::DynRankView< worksetJacobianValueType, worksetJacobianProperties...> worksetJacobians, const ordinal_type worksetFaceOrd, const shards::CellTopology parentCell) |
Computes non-normalized tangent vector pairs to physical faces in a face workset ; (see Subcell worksets for definition of face worksets). More... | |
template<typename faceTanValueType , class... faceTanProperties, typename worksetJacobianValueType , class... worksetJacobianProperties, typename faceOrdValueType , class... faceOrdProperties> | |
static void | getPhysicalFaceTangents (Kokkos::DynRankView< faceTanValueType, faceTanProperties...> faceTanU, Kokkos::DynRankView< faceTanValueType, faceTanProperties...> faceTanV, const Kokkos::DynRankView< worksetJacobianValueType, worksetJacobianProperties...> worksetJacobians, const Kokkos::DynRankView< faceOrdValueType, faceOrdProperties...> worksetFaceOrds, const shards::CellTopology parentCell) |
Computes non-normalized tangent vector pairs to physical faces in a face workset ; (see Subcell worksets for definition of face worksets). More... | |
template<typename sideNormalValueType , class... sideNormalProperties, typename worksetJacobianValueType , class... worksetJacobianProperties> | |
static void | getPhysicalSideNormals (Kokkos::DynRankView< sideNormalValueType, sideNormalProperties...> sideNormals, const Kokkos::DynRankView< worksetJacobianValueType, worksetJacobianProperties...> worksetJacobians, const ordinal_type worksetSideOrd, const shards::CellTopology parentCell) |
Computes non-normalized normal vectors to physical sides in a side workset . More... | |
template<typename sideNormalValueType , class... sideNormalProperties, typename worksetJacobianValueType , class... worksetJacobianProperties, typename edgeOrdValueType , class... edgeOrdProperties> | |
static void | getPhysicalSideNormals (Kokkos::DynRankView< sideNormalValueType, sideNormalProperties...> sideNormals, const Kokkos::DynRankView< worksetJacobianValueType, worksetJacobianProperties...> worksetJacobians, const Kokkos::DynRankView< edgeOrdValueType, edgeOrdProperties...> worksetSideOrds, const shards::CellTopology parentCell) |
Computes non-normalized normal vectors to physical sides in a side workset . It is similar to the CellTools::getPhysicalSideNormals function above, with the difference that the side ordinal can change from point to point, and it is provided by the rank-2 input array worksetSideOrds, indexed by (C,P). More... | |
template<typename faceNormalValueType , class... faceNormalProperties, typename worksetJacobianValueType , class... worksetJacobianProperties> | |
static void | getPhysicalFaceNormals (Kokkos::DynRankView< faceNormalValueType, faceNormalProperties...> faceNormals, const Kokkos::DynRankView< worksetJacobianValueType, worksetJacobianProperties...> worksetJacobians, const ordinal_type worksetFaceOrd, const shards::CellTopology parentCell) |
Computes non-normalized normal vectors to physical faces in a face workset ; (see Subcell worksets for definition of face worksets). More... | |
template<typename faceNormalValueType , class... faceNormalProperties, typename worksetJacobianValueType , class... worksetJacobianProperties, typename faceOrdValueType , class... faceOrdProperties> | |
static void | getPhysicalFaceNormals (Kokkos::DynRankView< faceNormalValueType, faceNormalProperties...> faceNormals, const Kokkos::DynRankView< worksetJacobianValueType, worksetJacobianProperties...> worksetJacobians, const Kokkos::DynRankView< faceOrdValueType, faceOrdProperties...> worksetFaceOrds, const shards::CellTopology parentCell) |
Computes non-normalized normal vectors to physical faces in a face workset ; (see Subcell worksets for definition of face worksets). It is similar to the CellTools::getPhysicalSideNormals function above, with the difference that the side ordinal can change from point to point, and it is provided by the rank-2 input array worksetSideOrds, indexed by (C,P). More... | |
template<typename PhysPointValueType , typename RefPointValueType , typename WorksetType , typename HGradBasisPtrType > | |
static void | mapToPhysicalFrame (PhysPointValueType physPoints, const RefPointValueType refPoints, const WorksetType worksetCell, const HGradBasisPtrType basis) |
Computes F, the reference-to-physical frame map. More... | |
template<typename PhysPointViewType , typename RefPointViewType , typename WorksetCellViewType > | |
static void | mapToPhysicalFrame (PhysPointViewType physPoints, const RefPointViewType refPoints, const WorksetCellViewType worksetCell, const shards::CellTopology cellTopo) |
Computes F, the reference-to-physical frame map. More... | |
template<typename refSubcellViewType , typename paramPointViewType > | |
static void | mapToReferenceSubcell (refSubcellViewType refSubcellPoints, const paramPointViewType paramPoints, const ordinal_type subcellDim, const ordinal_type subcellOrd, const shards::CellTopology parentCell) |
Computes parameterization maps of 1- and 2-subcells of reference cells. More... | |
template<typename refSubcellViewType , typename paramPointViewType > | |
static void | mapToReferenceSubcell (refSubcellViewType refSubcellPoints, const paramPointViewType paramPoints, const typename RefSubcellParametrization< DeviceType >::ConstViewType subcellParametrization, const ordinal_type subcellOrd) |
Computes parameterization maps of 1- and 2-subcells of reference cells. More... | |
template<typename refSubcellViewType , typename paramPointViewType , typename ordViewType > | |
static void | mapToReferenceSubcell (refSubcellViewType refSubcellPoints, const paramPointViewType paramPoints, const typename RefSubcellParametrization< DeviceType >::ConstViewType subcellParametrization, const ordViewType subcellOrd) |
Computes parameterization maps of 1- and 2-subcells of reference cells. More... | |
template<typename refPointValueType , class... refPointProperties, typename physPointValueType , class... physPointProperties, typename worksetCellValueType , class... worksetCellProperties> | |
static void | mapToReferenceFrame (Kokkos::DynRankView< refPointValueType, refPointProperties...> refPoints, const Kokkos::DynRankView< physPointValueType, physPointProperties...> physPoints, const Kokkos::DynRankView< worksetCellValueType, worksetCellProperties...> worksetCell, const shards::CellTopology cellTopo) |
Computes , the inverse of the reference-to-physical frame map using a default initial guess. More... | |
template<typename refPointValueType , class... refPointProperties, typename initGuessValueType , class... initGuessProperties, typename physPointValueType , class... physPointProperties, typename worksetCellValueType , class... worksetCellProperties, typename HGradBasisPtrType > | |
static void | mapToReferenceFrameInitGuess (Kokkos::DynRankView< refPointValueType, refPointProperties...> refPoints, const Kokkos::DynRankView< initGuessValueType, initGuessProperties...> initGuess, const Kokkos::DynRankView< physPointValueType, physPointProperties...> physPoints, const Kokkos::DynRankView< worksetCellValueType, worksetCellProperties...> worksetCell, const HGradBasisPtrType basis) |
Computation of , the inverse of the reference-to-physical frame map using user-supplied initial guess. More... | |
template<typename refPointValueType , class... refPointProperties, typename initGuessValueType , class... initGuessProperties, typename physPointValueType , class... physPointProperties, typename worksetCellValueType , class... worksetCellProperties> | |
static void | mapToReferenceFrameInitGuess (Kokkos::DynRankView< refPointValueType, refPointProperties...> refPoints, const Kokkos::DynRankView< initGuessValueType, initGuessProperties...> initGuess, const Kokkos::DynRankView< physPointValueType, physPointProperties...> physPoints, const Kokkos::DynRankView< worksetCellValueType, worksetCellProperties...> worksetCell, const shards::CellTopology cellTopo) |
Computation of , the inverse of the reference-to-physical frame map using user-supplied initial guess. More... | |
template<typename subcvCoordValueType , class... subcvCoordProperties, typename cellCoordValueType , class... cellCoordProperties> | |
static void | getSubcvCoords (Kokkos::DynRankView< subcvCoordValueType, subcvCoordProperties...> subcvCoords, const Kokkos::DynRankView< cellCoordValueType, cellCoordProperties...> cellCoords, const shards::CellTopology primaryCell) |
Computes coordinates of sub-control volumes in each primary cell. More... | |
template<typename PointViewType > | |
static bool | checkPointInclusion (const PointViewType point, const shards::CellTopology cellTopo, const typename ScalarTraits< typename PointViewType::value_type >::scalar_type thres=threshold< typename ScalarTraits< typename PointViewType::value_type >::scalar_type >()) |
Checks if a point belongs to a reference cell. More... | |
template<unsigned cellTopologyKey, typename OutputViewType , typename InputViewType > | |
static void | checkPointwiseInclusion (OutputViewType inCell, const InputViewType points, const typename ScalarTraits< typename InputViewType::value_type >::scalar_type thresh=threshold< typename ScalarTraits< typename InputViewType::value_type >::scalar_type >()) |
Checks every point for inclusion in the reference cell of a given topology. The points can belong to a global set and stored in a rank-2 view (P,D) , or to multiple sets indexed by a cell ordinal and stored in a rank-3 view (C,P,D). The cell topology key is a template argument. Requires cell topology with a reference cell. More... | |
template<typename InCellViewType , typename PointViewType > | |
static void | checkPointwiseInclusion (InCellViewType inCell, const PointViewType points, const shards::CellTopology cellTopo, const typename ScalarTraits< typename PointViewType::value_type >::scalar_type thres=threshold< typename ScalarTraits< typename PointViewType::value_type >::scalar_type >()) |
Checks every point in multiple sets indexed by a cell ordinal for inclusion in the reference cell of a given topology. Requires cell topology with a reference cell. More... | |
template<typename inCellValueType , class... inCellProperties, typename pointValueType , class... pointProperties, typename cellWorksetValueType , class... cellWorksetProperties> | |
static void | checkPointwiseInclusion (Kokkos::DynRankView< inCellValueType, inCellProperties...> inCell, const Kokkos::DynRankView< pointValueType, pointProperties...> points, const Kokkos::DynRankView< cellWorksetValueType, cellWorksetProperties...> cellWorkset, const shards::CellTopology cellTopo, const typename ScalarTraits< pointValueType >::scalar_type thres=threshold< typename ScalarTraits< pointValueType >::scalar_type >()) |
Checks every points for inclusion in physical cells from a cell workset. The points can belong to a global set and stored in a rank-2 (P,D) view, or to multiple sets indexed by a cell ordinal and stored in a rank-3 (C,P,D) view cells in a cell workset. More... | |
Private Types | |
using | ExecSpaceType = typename DeviceType::execution_space |
using | MemSpaceType = typename DeviceType::memory_space |
Static Private Member Functions | |
template<typename outputValueType , typename pointValueType > | |
static Teuchos::RCP< Basis < DeviceType, outputValueType, pointValueType > > | createHGradBasis (const shards::CellTopology cellTopo) |
Generates default HGrad basis based on cell topology. More... | |
A stateless class for operations on cell data. Provides methods for:
Definition at line 77 of file Intrepid2_CellTools.hpp.
|
static |
Allocates and returns a Data container suitable for storing determinants corresponding to the Jacobians in the Data container provided.
jacobian | [in] - data with shape (C,P,D,D), as returned by CellGeometry::allocateJacobianData() |
Definition at line 92 of file Intrepid2_CellToolsDefJacobian.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::getExtents(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView1(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView3(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView4(), and Intrepid2::Data< DataScalar, DeviceType >::getVariationTypes().
|
static |
Allocates and returns a Data container suitable for storing inverses corresponding to the Jacobians in the Data container provided.
jacobian | [in] - data with shape (C,P,D,D), as returned by CellGeometry::allocateJacobianData() |
Definition at line 126 of file Intrepid2_CellToolsDefJacobian.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::getExtents(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView1(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView2(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView3(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView4(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingViewRank(), Intrepid2::Data< DataScalar, DeviceType >::getVariationTypes(), and INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE.
|
static |
Checks if a point belongs to a reference cell.
Requires cell topology with a reference cell.
point | [in] - rank-1 view (D) of the point tested for inclusion |
cellTopo | [in] - cell topology |
threshold | [in] - "tightness" of the inclusion test |
Definition at line 37 of file Intrepid2_CellToolsDefInclusion.hpp.
|
static |
Checks every point for inclusion in the reference cell of a given topology. The points can belong to a global set and stored in a rank-2 view (P,D) , or to multiple sets indexed by a cell ordinal and stored in a rank-3 view (C,P,D). The cell topology key is a template argument. Requires cell topology with a reference cell.
inCell | [out] - rank-1 view (P) or rank-2 view (C,P). On return, its entries will be set to 1 or 0 depending on whether points are included in cells |
point | [in] - rank-2 view (P,D) or rank-3 view (C,P,D) with reference coordinates of the points tested for inclusion |
threshold | [in] - "tightness" of the inclusion test |
Definition at line 131 of file Intrepid2_CellToolsDefInclusion.hpp.
|
static |
Checks every point in multiple sets indexed by a cell ordinal for inclusion in the reference cell of a given topology. Requires cell topology with a reference cell.
inRefCell | [out] - rank-2 view (C,P) with results from the pointwise inclusion test |
refPoints | [in] - rank-3 view (C,P,D) |
cellTopo | [in] - cell topology |
threshold | [in] - "tightness" of the inclusion test |
|
static |
Checks every points for inclusion in physical cells from a cell workset. The points can belong to a global set and stored in a rank-2 (P,D) view, or to multiple sets indexed by a cell ordinal and stored in a rank-3 (C,P,D) view cells in a cell workset.
inCell | [out] - rank-2 view (P,D) with results from the pointwise inclusion test |
points | [in] - rank-2 view (P,D) or rank-3 view (C,P,D) with the physical points |
cellWorkset | [in] - rank-3 view with dimensions (C,N,D) with the nodes of the cell workset |
cellTopo | [in] - cell topology |
threshold | [in] - tolerance for inclusion tests on the input points |
Definition at line 216 of file Intrepid2_CellToolsDefInclusion.hpp.
References Intrepid2::RealSpaceTools< DeviceType >::clone().
|
inlinestaticprivate |
Generates default HGrad basis based on cell topology.
cellTopo | [in] - cell topology |
Definition at line 100 of file Intrepid2_CellTools.hpp.
|
static |
Computes non-normalized tangent vectors to physical edges in an edge workset ; (see Subcell worksets for definition of edge worksets).
For every edge in the workset the tangents are computed at the points that are images of points from R=[-1,1] on edge . Returns rank-3 array with dimensions (C,P,D1), D1=2 or D1=3 such that
In this formula:
worksetJacobians
must contain the values of , where , i.e., Jacobians of the parent cells evaluated at points that are located on reference edge having the same local ordinal as the edges in the workset.edgeTangents | [out] - rank-3 array (C,P,D1) with tangents on workset edges |
worksetJacobians | [in] - rank-4 array (C,P,D1,D1) with Jacobians evaluated at ref. edge points |
worksetEdgeOrd | [in] - edge ordinal, relative to ref. cell, of the edge workset |
parentCell | [in] - cell topology of the parent reference cell |
Definition at line 372 of file Intrepid2_CellToolsDefNodeInfo.hpp.
References Intrepid2::RealSpaceTools< DeviceType >::matvec().
Referenced by Intrepid2::FunctionSpaceTools< DeviceType >::computeEdgeMeasure().
|
static |
Computes non-normalized tangent vectors to physical edges in an edge workset ; (see Subcell worksets for definition of edge worksets).
It is similar to the CellTools::getPhysicalEdgeTangents function above, with the difference that the edge ordinal can change from point to point, and it is provided by the rank-2 input array worksetEdgeOrds, indexed by (C,P).
edgeTangents | [out] - rank-3 array (C,P,D1) with tangents on workset edges |
worksetJacobians | [in] - rank-4 array (C,P,D1,D1) with Jacobians evaluated at ref. edge points |
worksetEdgeOrds | [in] - rank-2 array (C,P) with edge ordinals, relative to ref. cell, of the edge workset |
parentCell | [in] - cell topology of the parent reference cell |
Definition at line 454 of file Intrepid2_CellToolsDefNodeInfo.hpp.
References Intrepid2::RefSubcellParametrization< DeviceType >::get(), and Intrepid2::RealSpaceTools< DeviceType >::matvec().
|
static |
Computes non-normalized normal vectors to physical faces in a face workset ; (see Subcell worksets for definition of face worksets).
For every face in the workset the normals are computed at the points that are images of points from the parametrization domain R on face . Returns rank-3 array with dimensions (C,P,D), D=3, such that
In this formula:
worksetJacobians
must contain the values of , where , i.e., Jacobians of the parent cells evaluated at points that are located on reference face having the same local ordinal as the faces in the workset.faceNormals | [out] - rank-3 array (C,P,D), normals at workset faces |
worksetJacobians | [in] - rank-4 array (C,P,D,D) with Jacobians at ref. face points |
worksetFaceOrd | [in] - face ordinal, relative to ref. cell, of the face workset |
parentCell | [in] - cell topology of the parent reference cell |
Definition at line 813 of file Intrepid2_CellToolsDefNodeInfo.hpp.
References Intrepid2::RealSpaceTools< DeviceType >::vecprod().
Referenced by Intrepid2::FunctionSpaceTools< DeviceType >::computeFaceMeasure().
|
static |
Computes non-normalized normal vectors to physical faces in a face workset ; (see Subcell worksets for definition of face worksets). It is similar to the CellTools::getPhysicalSideNormals function above, with the difference that the side ordinal can change from point to point, and it is provided by the rank-2 input array worksetSideOrds, indexed by (C,P).
faceNormals | [out] - rank-3 array (C,P,D), normals at workset faces |
worksetJacobians | [in] - rank-4 array (C,P,D,D) with Jacobians at ref. face points |
worksetFaceOrds | [in] - rank-2 array (C,P) with face ordinals, relative to ref. cell, of the face workset |
parentCell | [in] - cell topology of the parent reference cell |
Definition at line 875 of file Intrepid2_CellToolsDefNodeInfo.hpp.
References Intrepid2::RealSpaceTools< DeviceType >::vecprod().
|
static |
Computes non-normalized tangent vector pairs to physical faces in a face workset ; (see Subcell worksets for definition of face worksets).
For every face in the workset the tangents are computed at the points that are images of points from the parametrization domain R on face . Returns 2 rank-3 arrays with dimensions (C,P,D), D=3 such that
In this formula:
worksetJacobians
must contain the values of , where , i.e., Jacobians of the parent cells evaluated at points that are located on reference face having the same local ordinal as the faces in the workset.faceTanU | [out] - rank-3 array (C,P,D), image of ref. face u-tangent at workset faces |
faceTanV | [out] - rank-3 array (C,P,D), image of ref. face u-tangent at workset faces |
worksetJacobians | [in] - rank-4 array (C,P,D,D) with Jacobians at ref. face points |
worksetFaceOrd | [in] - face ordinal, relative to ref. cell, of the face workset |
parentCell | [in] - cell topology of the parent reference cell |
Definition at line 518 of file Intrepid2_CellToolsDefNodeInfo.hpp.
References Intrepid2::RealSpaceTools< DeviceType >::matvec().
|
static |
Computes non-normalized tangent vector pairs to physical faces in a face workset ; (see Subcell worksets for definition of face worksets).
It is similar to the CellTools::getPhysicalFaceTangents function above, with the difference that the face ordinal can change from point to point, and it is provided by the rank-2 input array worksetFaceOrds, indexed by (C,P).
faceTanU | [out] - rank-3 array (C,P,D), image of ref. face u-tangent at workset faces |
faceTanV | [out] - rank-3 array (C,P,D), image of ref. face u-tangent at workset faces |
worksetJacobians | [in] - rank-4 array (C,P,D,D) with Jacobians at ref. face points |
worksetFaceOrds | [in] - rank-2 array (C,P) with face ordinals, relative to ref. cell, of the face workset |
parentCell | [in] - cell topology of the parent reference cell |
Definition at line 614 of file Intrepid2_CellToolsDefNodeInfo.hpp.
References Intrepid2::RefSubcellParametrization< DeviceType >::get(), and Intrepid2::RealSpaceTools< DeviceType >::matvec().
|
static |
Computes non-normalized normal vectors to physical sides in a side workset .
For every side in the workset the normals are computed at the points that are images of points from the parametrization domain R on side . A side is defined as a subcell of dimension one less than that of its parent cell. Therefore, sides of 2D cells are 1-subcells (edges) and sides of 3D cells are 2-subcells (faces).
Returns rank-3 array with dimensions (C,P,D), D = 2 or 3, such that
In this formula:
worksetJacobians
must contain the values of , where , i.e., Jacobians of the parent cells evaluated at points that are located on reference side having the same local ordinal as the sides in the workset.sideNormals | [out] - rank-3 array (C,P,D), normals at workset sides |
worksetJacobians | [in] - rank-4 array (C,P,D,D) with Jacobians at ref. side points |
worksetSideOrd | [in] - side ordinal, relative to ref. cell, of the side workset |
parentCell | [in] - cell topology of the parent reference cell |
Definition at line 715 of file Intrepid2_CellToolsDefNodeInfo.hpp.
Referenced by Intrepid2::CubatureControlVolumeSide< DeviceType, pointValueType, weightValueType >::getCubature().
|
static |
Computes non-normalized normal vectors to physical sides in a side workset . It is similar to the CellTools::getPhysicalSideNormals function above, with the difference that the side ordinal can change from point to point, and it is provided by the rank-2 input array worksetSideOrds, indexed by (C,P).
sideNormals | [out] - rank-3 array (C,P,D), normals at workset sides |
worksetJacobians | [in] - rank-4 array (C,P,D,D) with Jacobians at ref. side points |
worksetSideOrds | [in] - rank-2 array (C,P) with side ordinals, relative to ref. cell, of the side workset |
parentCell | [in] - cell topology of the parent reference cell |
Definition at line 766 of file Intrepid2_CellToolsDefNodeInfo.hpp.
|
static |
Computes the Cartesian coordinates of reference cell barycenter.
Requires cell topology with a reference cell.
cellCenter | [out] - coordinates of the specified reference cell center |
cell | [in] - cell topology |
Definition at line 36 of file Intrepid2_CellToolsDefNodeInfo.hpp.
References Intrepid2::RefCellCenter< DeviceType >::get().
Referenced by Intrepid2::Basis_HVOL_C0_FEM< DeviceType, outputValueType, pointValueType >::Basis_HVOL_C0_FEM().
|
static |
Computes constant tangent vectors to edges of 2D or 3D reference cells.
Returns rank-1 array with dimension (D), D=2 or D=3; such that
where is the parametrization map of the specified reference edge , given by
The length of computed edge tangents is one-half the length of their associated edges:
Because the edges of all reference cells are always affine images of [-1,1], the edge tangent is constant vector field.
refEdgeTangent | [out] - rank-1 array (D) with the edge tangent; D = cell dimension |
edgeOrd | [in] - ordinal of the edge whose tangent is computed |
parentCell | [in] - cell topology of the parent reference cell |
Definition at line 211 of file Intrepid2_CellToolsDefNodeInfo.hpp.
References Intrepid2::RefSubcellParametrization< DeviceType >::get().
Referenced by Intrepid2::Basis_HCURL_TET_In_FEM< DeviceType, outputValueType, pointValueType >::Basis_HCURL_TET_In_FEM(), Intrepid2::Basis_HCURL_TRI_In_FEM< DeviceType, outputValueType, pointValueType >::Basis_HCURL_TRI_In_FEM(), Intrepid2::ProjectionTools< DeviceType >::getHCurlBasisCoeffs(), Intrepid2::ProjectionTools< DeviceType >::getHGradBasisCoeffs(), and Intrepid2::ProjectionTools< DeviceType >::getL2BasisCoeffs().
|
static |
Computes constant normal vectors to faces of 3D reference cell.
Returns rank-1 array with dimension (D), D=3 such that
where is the parametrization map of the specified reference face given by
and
The length of computed face normals is proportional to face area:
Because the faces of all reference cells are always affine images of R , the coordinate functions of the parametrization map are linear and the face normal is a constant vector.
refFaceNormal | [out] - rank-1 array (D) with (constant) face normal |
faceOrd | [in] - ordinal of the face whose normal is computed |
parentCell | [in] - cell topology of the parent reference cell |
Definition at line 335 of file Intrepid2_CellToolsDefNodeInfo.hpp.
References Intrepid2::RealSpaceTools< DeviceType >::vecprod().
Referenced by Intrepid2::ProjectionTools< DeviceType >::getHCurlBasisCoeffs(), and Intrepid2::ProjectionTools< DeviceType >::getHGradBasisCoeffs().
|
static |
Computes pairs of constant tangent vectors to faces of a 3D reference cells.
Returns 2 rank-1 arrays with dimension (D), D=3, such that
where is the parametrization map of the specified reference face given by
and
Because the faces of all reference cells are always affine images of R , the coordinate functions of the parametrization map are linear and the face tangents are constant vectors.
refFaceTanU | [out] - rank-1 array (D) with (constant) tangent in u-direction |
refFaceTanV | [out] - rank-1 array (D) with (constant) tangent in v-direction |
faceOrd | [in] - ordinal of the face whose tangents are computed |
parentCell | [in] - cell topology of the parent 3D reference cell |
Definition at line 249 of file Intrepid2_CellToolsDefNodeInfo.hpp.
References Intrepid2::RefSubcellParametrization< DeviceType >::get().
Referenced by Intrepid2::Basis_HCURL_TET_In_FEM< DeviceType, outputValueType, pointValueType >::Basis_HCURL_TET_In_FEM().
|
static |
Retrieves the Cartesian coordinates of a reference cell node.
Returns Cartesian coordinates of a reference cell node. Requires cell topology with a reference cell. Node coordinates are always returned as an (x,y,z)-triple regardlesss of the actual topological cell dimension. The unused coordinates are set to zero, e.g., node 0 of Line<2> is returned as {-1,0,0}.
cellNode | [out] - coordinates of the specified reference vertex |
cell | [in] - cell topology of the cell |
vertexOrd | [in] - node ordinal |
Definition at line 139 of file Intrepid2_CellToolsDefNodeInfo.hpp.
References Intrepid2::RefCellNodes< DeviceType >::get().
|
static |
Computes constant normal vectors to sides of 2D or 3D reference cells.
A side is defined as a subcell of dimension one less than that of its parent cell. Therefore, sides of 2D cells are 1-subcells (edges) and sides of 3D cells are 2-subcells (faces).
Returns rank-1 array with dimension (D), D = 2 or 3 such that
where , and is the parametrization map of the specified reference side given by
For sides of 2D cells R=[-1,1] and for sides of 3D cells
For 3D cells the length of computed side normals is proportional to side area:
For 2D cells the length of computed side normals is proportional to side length:
Because the sides of all reference cells are always affine images of R , the coordinate functions of the parametrization maps are linear and the side normal is a constant vector.
refSideNormal | [out] - rank-1 array (D) with (constant) side normal |
sideOrd | [in] - ordinal of the side whose normal is computed |
parentCell | [in] - cell topology of the parent reference cell |
Definition at line 294 of file Intrepid2_CellToolsDefNodeInfo.hpp.
Referenced by Intrepid2::Basis_HDIV_TET_In_FEM< DeviceType, outputValueType, pointValueType >::Basis_HDIV_TET_In_FEM(), Intrepid2::Basis_HDIV_TRI_In_FEM< DeviceType, outputValueType, pointValueType >::Basis_HDIV_TRI_In_FEM(), Intrepid2::ProjectionTools< DeviceType >::getHDivBasisCoeffs(), and Intrepid2::ProjectionTools< DeviceType >::getL2BasisCoeffs().
|
static |
Retrieves the Cartesian coordinates of all nodes of a reference subcell.
Returns rank-2 array with the Cartesian coordinates of the nodes of the specified reference cell subcell. Requires cell topology with a reference cell.
subcellNodes | [out] - rank-2 (N,D) array with the Cartesian coordinates of the reference subcell nodes |
subcellDim | [in] - dimension of the subcell; 0 <= subcellDim <= parentCell dimension |
subcellOrd | [in] - ordinal of the subcell |
parentCell | [in] - topology of the cell that owns the subcell |
subcellDim
= dimension of the parentCell
this method returns the Cartesian coordinates of the nodes of the reference cell itself. Note that this requires subcellOrd=0
. Definition at line 169 of file Intrepid2_CellToolsDefNodeInfo.hpp.
|
static |
Retrieves the Cartesian coordinates of all vertices of a reference subcell.
Returns rank-2 array with the Cartesian coordinates of the vertices of the specified reference cell subcell. Requires cell topology with a reference cell.
subcellVertices | [out] - rank-2 (V,D) array with the Cartesian coordinates of the reference subcell vertices |
subcellDim | [in] - dimension of the subcell; 0 <= subcellDim <= parentCell dimension |
subcellOrd | [in] - ordinal of the subcell |
parentCell | [in] - topology of the cell that owns the subcell |
subcellDim
= dimension of the parentCell
this method returns the Cartesian coordinates of the vertices of the reference cell itself. Note that this requires subcellOrd=0. Definition at line 103 of file Intrepid2_CellToolsDefNodeInfo.hpp.
|
static |
Retrieves the Cartesian coordinates of a reference cell vertex.
Requires cell topology with a reference cell.
cellVertex | [out] - coordinates of the specified reference cell vertex |
cell | [in] - cell topology of the cell |
vertexOrd | [in] - vertex ordinal |
Definition at line 68 of file Intrepid2_CellToolsDefNodeInfo.hpp.
References Intrepid2::RefCellNodes< DeviceType >::get().
Referenced by Intrepid2::ProjectionStruct< DeviceType, ValueType >::createHGradProjectionStruct(), and Intrepid2::ProjectionStruct< DeviceType, ValueType >::createL2ProjectionStruct().
|
static |
Computes coordinates of sub-control volumes in each primary cell.
To build the system of equations for the control volume finite element method we need to compute geometric data for integration over control volumes. A control volume is polygon or polyhedron that surrounds a primary cell node and has vertices that include the surrounding primary cells' barycenter, edge midpoints, and face midpoints if in 3-d.
When using element-based assembly of the discrete equations over the primary mesh, a single element will contain a piece of each control volume surrounding each of the primary cell nodes. This piece of control volume (sub-control volume) is always a quadrilateral in 2-d and a hexahedron in 3-d.
In 2-d the sub-control volumes are defined in the following way:
Quadrilateral primary element: O________M________O | | | | 3 | 2 | B = cell barycenter | | | O = primary cell nodes M________B________M M = cell edge midpoints | | | | 0 | 1 | sub-control volumes 0, 1, 2, 3 | | | O________M________O Triangle primary element: O / \ / \ B = cell barycenter / \ O = primary cell nodes M 2 M M = cell edge midpoints / \ / \ / \ B / \ sub-control volumes 0, 1, 2 / | \ / 0 | 1 \ O________M________O
In 3-d the sub-control volumes are defined by the primary cell face centers and edge midpoints. The eight sub-control volumes for a hexahedron are shown below:
O__________E__________O /| /| /| E_|________F_|________E | /| | /| | /| | O_|_|______E_|_|______O | | O = primary cell nodes | | E------|-|-F------|-|-E B = cell barycenter | |/| | |/| | |/| F = cell face centers | F-|------|-B-|------|-F | E = cell edge midpoints |/| | |/| | |/| | E_|_|______F_|_|______E | | | | O------|-|-E------|-|-O | |/ | |/ | |/ | E--------|-F--------|-E |/ |/ |/ O__________E__________O
subCVCoords | [out] - array containing sub-control volume coordinates |
cellCoords | [in] - array containing coordinates of primary cells |
primaryCell | [in] - primary cell topology |
Definition at line 340 of file Intrepid2_CellToolsDefControlVolume.hpp.
Referenced by Intrepid2::CubatureControlVolumeBoundary< DeviceType, pointValueType, weightValueType >::getCubature(), Intrepid2::CubatureControlVolume< DeviceType, pointValueType, weightValueType >::getCubature(), and Intrepid2::CubatureControlVolumeSide< DeviceType, pointValueType, weightValueType >::getCubature().
|
inlinestatic |
Checks if a cell topology has reference cell.
cell | [in] - cell topology |
Definition at line 88 of file Intrepid2_CellTools.hpp.
References Intrepid2::RefSubcellParametrization< DeviceType >::isSupported().
|
static |
Computes F, the reference-to-physical frame map.
There are 2 use cases:
For a single point set in a rank-2 array (P,D) returns a rank-3 (C,P,D) array such that
For multiple point sets in a rank-3 (C,P,D) array returns a rank-3 (C,P,D) array such that
This corresponds to mapping multiple sets of reference points to a matching number of physical cells.
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.
refPoints
represents an arbitrary set of points in the reference frame that are not required to be in the reference cell corresponding to the specified cell topology. As a result, the images of these points under a given reference-to-physical map are not necessarily contained in the physical cell that is the image of the reference cell under that map. CellTools provides several inclusion tests methods to check whether or not the points are inside a reference cell.physPoints | [out] - rank-3 array with dimensions (C,P,D) with the images of the ref. points |
refPoints | [in] - rank-3/2 array with dimensions (C,P,D)/(P,D) with points in reference frame |
cellWorkset | [in] - rank-3 container with logical dimensions (C,N,D) with the nodes of the cell workset |
basis | [in] - pointer to HGrad basis used in reference-to-physical cell mapping |
Referenced by Intrepid2::CellTools< DeviceType >::mapToPhysicalFrame(), and Intrepid2::ProjectedGeometry< spaceDim, PointScalar, DeviceType >::projectOntoHGRADBasis().
|
inlinestatic |
Computes F, the reference-to-physical frame map.
There are 2 use cases:
For a single point set in a rank-2 array (P,D) returns a rank-3 (C,P,D) array such that
For multiple point sets in a rank-3 (C,P,D) array returns a rank-3 (C,P,D) array such that
This corresponds to mapping multiple sets of reference points to a matching number of physical cells.
Requires cell topology with a reference cell. See Section Reference-to-physical cell mapping for definition of the mapping function. Presently supported cell topologies are
Line<2>
Triangle<3>
, Triangle<6>
, Quadrilateral<4>
, Quadrilateral<9>
Tetrahedron<4>
, Tetrahedron<10>
, Hexahedron<8>
, Hexahedron<27>
refPoints
represents an arbitrary set of points in the reference frame that are not required to be in the reference cell corresponding to the specified cell topology. As a result, the images of these points under a given reference-to-physical map are not necessarily contained in the physical cell that is the image of the reference cell under that map. CellTools provides several inclusion tests methods to check whether or not the points are inside a reference cell.physPoints | [out] - rank-3 array with dimensions (C,P,D) with the images of the ref. points |
refPoints | [in] - rank-3/2 array with dimensions (C,P,D)/(P,D) with points in reference frame |
cellWorkset | [in] - rank-3 array with dimensions (C,N,D) with the nodes of the cell workset |
cellTopo | [in] - cell topology of the cells stored in cellWorkset |
Definition at line 1052 of file Intrepid2_CellTools.hpp.
References Intrepid2::CellTools< DeviceType >::mapToPhysicalFrame().
|
static |
Computes , the inverse of the reference-to-physical frame map using a default initial guess.
Applies for all cells in a cell workset to multiple point sets having the same number of points, indexed by cell ordinal, and stored in a rank-3 (C,P,D) array. Returns a rank-3 (C,P,D) array such that
Requires cell topology with a reference cell. See Section Reference-to-physical cell mapping for definition of the mapping function. Presently supported cell topologies are
Line<2>
Triangle<3>
, Triangle<6>
, Quadrilateral<4>
, Quadrilateral<9>
Tetrahedron<4>
, Tetrahedron<10>
, Hexahedron<8>
, Hexahedron<27>
Line
topologies: line center (0) Triangle
topologies: the point (1/3, 1/3) Quadrilateral
topologies: the point (0, 0) Tetrahedron
topologies: the point (1/6, 1/6, 1/6) Hexahedron
topologies: the point (0, 0, 0) Wedge
topologies: the point (1/2, 1/2, 0) For some cells with extended topologies, these initial guesses may not be good enough for Newton's method to converge in the allotted number of iterations. A version of this method with user-supplied initial guesses is also available.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.refPoints | [out] - rank-3 array with dimensions (C,P,D) with the images of the physical points |
physPoints | [in] - rank-3 array with dimensions (C,P,D) with points in physical frame |
worksetCell | [in] - rank-3 array with dimensions (C,N,D) with the nodes of the cell workset |
cellTopo | [in] - cell topology of the cells stored in cellWorkset |
Definition at line 39 of file Intrepid2_CellToolsDefPhysToRef.hpp.
Referenced by Intrepid2::PointTools::getWarpBlendLatticeTetrahedron(), and Intrepid2::PointTools::getWarpBlendLatticeTriangle().
|
static |
Computation of , the inverse of the reference-to-physical frame map using user-supplied initial guess.
Applies for all cells in a cell workset to multiple point sets having the same number of points, indexed by cell ordinal, and stored in a rank-3 (C,P,D) array. Returns a rank-3 (C,P,D) array such that
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.
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.refPoints | [out] - rank-3/2 array with dimensions (C,P,D)/(P,D) with the images of the physical points |
initGuess | [in] - rank-3/2 array with dimensions (C,P,D)/(P,D) with the initial guesses for each point |
physPoints | [in] - rank-3/2 array with dimensions (C,P,D)/(P,D) with points in physical frame |
worksetCell | [in] - rank-3 array with dimensions (C,N,D) with the nodes of the cell workset |
basis | [in] - pointer to HGrad basis used for reference to physical cell mapping |
Definition at line 90 of file Intrepid2_CellToolsDefPhysToRef.hpp.
References Intrepid2::Parameters::MaxNewton.
Referenced by Intrepid2::CellTools< DeviceType >::mapToReferenceFrameInitGuess().
|
inlinestatic |
Computation of , the inverse of the reference-to-physical frame map using user-supplied initial guess.
Applies for all cells in a cell workset to multiple point sets having the same number of points, indexed by cell ordinal, and stored in a rank-3 (C,P,D) array. Returns a rank-3 (C,P,D) array such that
Requires cell topology with a reference cell. See Section Reference-to-physical cell mapping for definition of the mapping function. Presently supported cell topologies are
Line<2>
Triangle<3>
, Triangle<6>
, Quadrilateral<4>
, Quadrilateral<9>
Tetrahedron<4>
, Tetrahedron<10>
, Hexahedron<8>
, Hexahedron<27>
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.refPoints | [out] - rank-3/2 array with dimensions (C,P,D)/(P,D) with the images of the physical points |
initGuess | [in] - rank-3/2 array with dimensions (C,P,D)/(P,D) with the initial guesses for each point |
physPoints | [in] - rank-3/2 array with dimensions (C,P,D)/(P,D) with points in physical frame |
worksetCell | [in] - rank-3 array with dimensions (C,N,D) with the nodes of the cell workset |
cellTopo | [in] - cell topology of the cells stored in cellWorkset |
Definition at line 1284 of file Intrepid2_CellTools.hpp.
References Intrepid2::CellTools< DeviceType >::mapToReferenceFrameInitGuess().
|
static |
Computes parameterization maps of 1- and 2-subcells of reference cells.
Applies , the parametrization map of a subcell from a given reference cell, to a set of points in the parametrization domain R of . Returns a rank-2 array with dimensions (P,D) where for 1-subcells:
for 2-subcells:
and
paramPoints
. Here c is ordinal of a parent cell, relative to subcell workset, and i is subcell ordinal, relative to a reference cell, of the subcell workset. See Section Subcell worksets for definition of subcell workset and Section Parametrization of physical 1- and 2-subcells for definition of parametrization maps.refSubcellPoints | [out] - rank-2 (P,D1) array with images of parameter space points |
paramPoints | [in] - rank-2 (P,D2) array with points in 1D or 2D parameter domain |
subcellDim | [in] - dimension of the subcell where points are mapped to |
subcellOrd | [in] - subcell ordinal |
parentCell | [in] - cell topology of the parent cell. |
Referenced by Intrepid2::Basis_HCURL_TET_In_FEM< DeviceType, outputValueType, pointValueType >::Basis_HCURL_TET_In_FEM(), Intrepid2::Basis_HCURL_TRI_In_FEM< DeviceType, outputValueType, pointValueType >::Basis_HCURL_TRI_In_FEM(), Intrepid2::Basis_HDIV_TET_In_FEM< DeviceType, outputValueType, pointValueType >::Basis_HDIV_TET_In_FEM(), Intrepid2::Basis_HDIV_TRI_In_FEM< DeviceType, outputValueType, pointValueType >::Basis_HDIV_TRI_In_FEM(), Intrepid2::Basis_HGRAD_TET_Cn_FEM< DeviceType, outputValueType, pointValueType >::Basis_HGRAD_TET_Cn_FEM(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::createHCurlProjectionStruct(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::createHDivProjectionStruct(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::createHGradProjectionStruct(), Intrepid2::ProjectionStruct< DeviceType, ValueType >::createL2ProjectionStruct(), Intrepid2::CubatureControlVolumeBoundary< DeviceType, pointValueType, weightValueType >::CubatureControlVolumeBoundary(), and Intrepid2::CubatureControlVolumeSide< DeviceType, pointValueType, weightValueType >::CubatureControlVolumeSide().
|
static |
Computes parameterization maps of 1- and 2-subcells of reference cells.
Overload of the previous function (see explanation above) where the subcell parametrization is used instead of passing the parent cell topology.
|
static |
Computes parameterization maps of 1- and 2-subcells of reference cells.
Overload of the previous function (see explanation above) where the subcellOrd is a rank-1 array (P).
|
static |
Computes the Jacobian matrix DF of the reference-to-physical frame map F.
There are two use cases:
For a single point set in a rank-2 array (P,D) returns a rank-4 (C,P,D,D) array such that
For multiple sets of reference points in a rank-3 (C,P,D) array returns rank-4 (C,P,D,D) array such that
Requires pointer to HGrad basis that defines reference to physical cell mapping. See Section \ref sec_cell_topology_ref_map_DF for definition of the Jacobian. \warning The points are not required to be in the reference cell associated with the specified cell topology. CellTools provides several inclusion tests methods to check whether or not the points are inside a reference cell.
jacobian | [out] - rank-4 array with dimensions (C,P,D,D) with the Jacobians |
points | [in] - rank-2/3 array with dimensions (P,D)/(C,P,D) with the evaluation points |
cellWorkset | [in] - rank-3 container with logical dimensions (C,N,D) with the nodes of the cell workset |
basis | [in] - HGrad basis for reference to physical cell mapping |
startCell | [in] - first cell index in cellWorkset for which we should compute the Jacobian; corresponds to the 0 index in Jacobian and/or points container. Default: 0. |
endCell | [in] - first cell index in cellWorkset that we do not process; endCell - startCell must equal the extent of the Jacobian container in dimension 0. Default: -1, a special value that indicates the extent of the cellWorkset should be used. |
Definition at line 799 of file Intrepid2_CellToolsDefJacobian.hpp.
Referenced by Intrepid2::CubatureControlVolumeBoundary< DeviceType, pointValueType, weightValueType >::getCubature(), Intrepid2::CubatureControlVolume< DeviceType, pointValueType, weightValueType >::getCubature(), Intrepid2::CubatureControlVolumeSide< DeviceType, pointValueType, weightValueType >::getCubature(), Intrepid2::CellTools< DeviceType >::setJacobian(), and Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::setJacobianDataPrivate().
|
static |
Computes the Jacobian matrix DF of the reference-to-physical frame map F.
There are two use cases:
For a single point set in a rank-2 array (P,D) returns a rank-4 (C,P,D,D) array such that
For multiple sets of reference points in a rank-3 (C,P,D) array returns rank-4 (C,P,D,D) array such that
Requires pointer to HGrad basis that defines reference to physical cell mapping. See Section \ref sec_cell_topology_ref_map_DF for definition of the Jacobian. \warning The points are not required to be in the reference cell associated with the specified cell topology. CellTools provides several inclusion tests methods to check whether or not the points are inside a reference cell.
jacobian | [out] - rank-4 array with dimensions (C,P,D,D) with the Jacobians |
cellWorkset | [in] - rank-3 container with logical dimensions (C,N,D) with the nodes of the cell workset |
gradients | [in] - rank-3/4 array with dimensions (N,P,D)/(C,N,P,D) with the gradients of the physical-to-cell mapping |
startCell | [in] - first cell index in cellWorkset for which we should compute the Jacobian; corresponds to the 0 index in Jacobian and/or points container. Default: 0. |
endCell | [in] - first cell index in cellWorkset that we do not process; endCell - startCell must equal the extent of the Jacobian container in dimension 0. Default: -1, a special value that indicates the extent of the cellWorkset should be used. |
Definition at line 771 of file Intrepid2_CellToolsDefJacobian.hpp.
|
inlinestatic |
Computes the Jacobian matrix DF of the reference-to-physical frame map F.
There are two use cases:
For a single point set in a rank-2 array (P,D) returns a rank-4 (C,P,D,D) array such that
For multiple sets of reference points in a rank-3 (C,P,D) array returns rank-4 (C,P,D,D) array such that
Requires cell topology with a reference cell. See Section \ref sec_cell_topology_ref_map_DF for definition of the Jacobian. \warning The points are not required to be in the reference cell associated with the specified cell topology. CellTools provides several inclusion tests methods to check whether or not the points are inside a reference cell.
jacobian | [out] - rank-4 array with dimensions (C,P,D,D) with the Jacobians |
points | [in] - rank-2/3 array with dimensions (P,D)/(C,P,D) with the evaluation points |
cellWorkset | [in] - rank-3 array with dimensions (C,N,D) with the nodes of the cell workset |
cellTopo | [in] - cell topology of the cells stored in cellWorkset |
Definition at line 284 of file Intrepid2_CellTools.hpp.
References Intrepid2::CellTools< DeviceType >::setJacobian().
|
static |
Computes the determinant of the Jacobian matrix DF of the reference-to-physical frame map F.
Returns rank-2 array with dimensions (C,P) such that
jacobianDet | [out] - rank-2 array with dimensions (C,P) with Jacobian determinants |
jacobian | [in] - rank-4 array with dimensions (C,P,D,D) with the Jacobians |
Definition at line 874 of file Intrepid2_CellToolsDefJacobian.hpp.
References Intrepid2::RealSpaceTools< DeviceType >::det().
Referenced by Intrepid2::CubatureControlVolumeBoundary< DeviceType, pointValueType, weightValueType >::getCubature(), and Intrepid2::CubatureControlVolume< DeviceType, pointValueType, weightValueType >::getCubature().
|
static |
Computes determinants corresponding to the Jacobians in the Data container provided.
jacobianDet | [out] - data with shape (C,P), as returned by CellTools::allocateJacobianDet() |
jacobian | [in] - data with shape (C,P,D,D), as returned by CellGeometry::allocateJacobianData() |
Definition at line 165 of file Intrepid2_CellToolsDefJacobian.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::blockPlusDiagonalLastNonDiagonal(), Intrepid2::RealSpaceTools< DeviceType >::det(), Intrepid2::Data< DataScalar, DeviceType >::extent_int(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView1(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView2(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView3(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView4(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingViewRank(), Intrepid2::Data< DataScalar, DeviceType >::getVariationTypes(), and INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE.
|
static |
Computes reciprocals of determinants corresponding to the Jacobians in the Data container provided.
jacobianDetInv | [out] - data with shape (C,P), as returned by CellTools::allocateJacobianDet() |
jacobian | [in] - data with shape (C,P,D,D), as returned by CellGeometry::allocateJacobianData() |
Definition at line 432 of file Intrepid2_CellToolsDefJacobian.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::allocateConstantData(), and Intrepid2::Data< DataScalar, DeviceType >::storeInPlaceQuotient().
|
static |
Multiplies the Jacobian with shape (C,P,D,D) by the reciprocals of the determinants, with shape (C,P), entrywise.
jacobianDividedByDet | [out] - data container with shape (C,P,D,D), as returned by CellTools::allocateJacobianInv() |
jacobian | [in] - data container with shape (C,P,D,D), as returned by CellTools::allocateJacobianInv() |
jacobianDetInv | [in] - data with shape (C,P,D,D), as returned by CellGeometry::allocateJacobianData() |
|
static |
Computes the inverse of the Jacobian matrix DF of the reference-to-physical frame map F.
Returns rank-4 array with dimensions (C,P,D,D) such that
jacobianInv | [out] - rank-4 array with dimensions (C,P,D,D) with the inverse Jacobians |
jacobian | [in] - rank-4 array with dimensions (C,P,D,D) with the Jacobians |
Definition at line 861 of file Intrepid2_CellToolsDefJacobian.hpp.
References Intrepid2::RealSpaceTools< DeviceType >::inverse().
|
static |
Computes determinants corresponding to the Jacobians in the Data container provided.
jacobianInv | [out] - data container with shape (C,P,D,D), as returned by CellTools::allocateJacobianInv() |
jacobian | [in] - data with shape (C,P,D,D), as returned by CellGeometry::allocateJacobianData() |
Definition at line 442 of file Intrepid2_CellToolsDefJacobian.hpp.
References Intrepid2::Data< DataScalar, DeviceType >::blockPlusDiagonalLastNonDiagonal(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView1(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView2(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView3(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView4(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingViewRank(), Intrepid2::Data< DataScalar, DeviceType >::getVariationTypes(), INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE, and Intrepid2::RealSpaceTools< DeviceType >::inverse().