Intrepid2
Public Types | Public Member Functions | Protected Types | Protected Attributes | List of all members
Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType > Class Template Reference

CellGeometry provides the nodes for a set of cells; has options that support efficient definition of uniform grids as well as options for arbitrary geometry, including curvilinear. More...

#include <Intrepid2_CellGeometry.hpp>

Public Types

enum  CellGeometryType {
  UNIFORM_GRID, TENSOR_GRID, EXTRUDED_GRID, FIRST_ORDER,
  HIGHER_ORDER
}
 
enum  SubdivisionStrategy {
  NO_SUBDIVISION, TWO_TRIANGLES_RIGHT, TWO_TRIANGLES_LEFT, FOUR_TRIANGLES,
  FIVE_TETRAHEDRA, SIX_TETRAHEDRA, SIX_PYRAMIDS
}
 
enum  HypercubeNodeOrdering { HYPERCUBE_NODE_ORDER_CLASSIC_SHARDS, HYPERCUBE_NODE_ORDER_TENSOR }
 

Public Member Functions

KOKKOS_INLINE_FUNCTION int numCellsPerGridCell (SubdivisionStrategy subdivisionStrategy) const
 Helper method that returns the number of cells into which each grid cell will be subdivided based on the input subdivision strategy. More...
 
Data< PointScalar, DeviceType > allocateJacobianDataPrivate (const ScalarView< PointScalar, DeviceType > &pointComponentView, const int &pointsPerCell, const int startCell, const int endCell) const
 Notionally-private method that provides a common interface for multiple public-facing allocateJacobianData() methods. (Marked as public due to compiler constraints.) More...
 
void setJacobianDataPrivate (Data< PointScalar, DeviceType > &jacobianData, const int &pointsPerCell, const Data< PointScalar, DeviceType > &refData, const int startCell, const int endCell) const
 Notionally-private method that provides a common interface for multiple public-facing setJacobianData() methods. (Marked as public due to compiler constraints.) More...
 
 CellGeometry (const Kokkos::Array< PointScalar, spaceDim > &origin, const Kokkos::Array< PointScalar, spaceDim > &domainExtents, const Kokkos::Array< int, spaceDim > &gridCellCounts, SubdivisionStrategy subdivisionStrategy=NO_SUBDIVISION, HypercubeNodeOrdering nodeOrdering=HYPERCUBE_NODE_ORDER_TENSOR)
 Uniform grid constructor, with optional subdivision into simplices. More...
 
 CellGeometry (const shards::CellTopology &cellTopo, ScalarView< int, DeviceType > cellToNodes, ScalarView< PointScalar, DeviceType > nodes, const bool claimAffine=false, const HypercubeNodeOrdering nodeOrdering=HYPERCUBE_NODE_ORDER_CLASSIC_SHARDS)
 Node-based constructor for straight-edged geometry. More...
 
 CellGeometry (Teuchos::RCP< Intrepid2::Basis< DeviceType, PointScalar, PointScalar > > basisForNodes, ScalarView< PointScalar, DeviceType > cellNodes)
 Node-based constructor for curvilinear geometry. More...
 
KOKKOS_INLINE_FUNCTION CellGeometry (const CellGeometry &cellGeometry)
 Copy constructor. More...
 
KOKKOS_INLINE_FUNCTION ~CellGeometry ()
 Destructor.
 
KOKKOS_INLINE_FUNCTION bool affine () const
 Returns true if Jacobian is constant within each cell.
 
TensorData< PointScalar,
DeviceType > 
allocateCellMeasure (const Data< PointScalar, DeviceType > &jacobianDet, const TensorData< PointScalar, DeviceType > &cubatureWeights) const
 Allocate a TensorData object appropriate for passing to computeCellMeasure(). More...
 
void computeCellMeasure (TensorData< PointScalar, DeviceType > &cellMeasure, const Data< PointScalar, DeviceType > &jacobianDet, const TensorData< PointScalar, DeviceType > &cubatureWeights) const
 Compute cell measures that correspond to provided Jacobian determinants and. More...
 
BasisPtr basisForNodes () const
 H^1 Basis used in the reference-to-physical transformation. Linear for straight-edged geometry; higher-order for curvilinear.
 
const shards::CellTopology & cellTopology () const
 The shards CellTopology for each cell within the CellGeometry object. Note that this is always a lowest-order CellTopology, even for higher-order curvilinear geometry. Higher-order geometry for CellGeometry is expressed in terms of the basis returned by basisForNodes().
 
KOKKOS_INLINE_FUNCTION ordinal_type cellToNode (ordinal_type cellIndex, ordinal_type cellNodeNumber) const
 Returns a global node number corresponding to the specified (cell, node) combination. If uniform grid (possibly subdivided), this number is computed dynamically; for more general meshes, this simply returns the result of a lookup in the cellToNodes container provided at construction.
 
KOKKOS_INLINE_FUNCTION
DataVariationType 
cellVariationType () const
 
KOKKOS_INLINE_FUNCTION int hypercubeComponentNodeNumber (int hypercubeNodeNumber, int d) const
 For hypercube vertex number hypercubeNodeNumber, returns the component node number in specified dimension d. This is either 0 or 1. (The mapping depends on the node ordering specified in nodeOrdering_.)
 
void initializeOrientations ()
 Initialize the internal orientations_ member with the orientations of each member cell. These are used for projected geometry, and have a logical shape (C).
 
KOKKOS_INLINE_FUNCTION size_t extent (const int &r) const
 Returns the logical extent of the container in the specified dimension; the shape of CellGeometry is always (C,N,D), where C is the number of cells, N is the number of nodes per cell (may be more than the number of vertices, in the case of curvilinear geometry), and D is the spatial dimension.
 
template<typename iType >
KOKKOS_INLINE_FUNCTION
std::enable_if
< std::is_integral< iType >
::value, int >::type 
extent_int (const iType &r) const
 Returns the logical extent of the container in the specified dimension as an int; the shape of CellGeometry is always (C,N,D), where C is the number of cells, N is the number of nodes per cell (may be more than the number of vertices, in the case of curvilinear geometry), and D is the spatial dimension.
 
KOKKOS_INLINE_FUNCTION
HypercubeNodeOrdering 
nodeOrderingForHypercubes () const
 Returns the node ordering used for hypercubes.
 
KOKKOS_INLINE_FUNCTION int numCells () const
 Returns the number of cells.
 
KOKKOS_INLINE_FUNCTION int numCellsInDimension (const int &dim) const
 For uniform grid and tensor grid CellGeometry, returns the number of cells in the specified component dimension. For other CellGeometry types, returns -1.
 
KOKKOS_INLINE_FUNCTION int numNodesPerCell () const
 Returns the number of nodes per cell; may be more than the number of vertices in the corresponding CellTopology, in the case of higher-order (curvilinear) geometry.
 
KOKKOS_INLINE_FUNCTION Orientation getOrientation (int &cellNumber) const
 Returns the orientation for the specified cell. Requires that initializeOrientations() has been called previously.
 
Data< Orientation, DeviceType > getOrientations ()
 Returns the orientations for all cells. Calls initializeOrientations() if it has not previously been called.
 
void orientations (ScalarView< Orientation, DeviceType > orientationsView, const int &startCell=0, const int &endCell=-1)
 Fills the provided container with the orientations for the specified cell range. Calls getOrientations() and copies the orientations from the Data container into the ScalarView container. More...
 
KOKKOS_INLINE_FUNCTION PointScalar gridCellCoordinate (const int &gridCellOrdinal, const int &localNodeNumber, const int &dim) const
 returns coordinate in dimension dim of the indicated node in the indicated grid cell
 
KOKKOS_INLINE_FUNCTION unsigned rank () const
 Returns the logical rank of this container. This is always 3.
 
KOKKOS_INLINE_FUNCTION int gridCellNodeForSubdivisionNode (const int &gridCellOrdinal, const int &subdivisionOrdinal, const int &subdivisionNodeNumber) const
 returns coordinate in dimension d for the indicated subdivision of the indicated grid cell
 
KOKKOS_INLINE_FUNCTION PointScalar subdivisionCoordinate (const int &gridCellOrdinal, const int &subdivisionOrdinal, const int &subdivisionNodeNumber, const int &d) const
 returns coordinate in dimension d for the indicated subdivision of the indicated grid cell
 
KOKKOS_INLINE_FUNCTION PointScalar operator() (const int &cell, const int &node, const int &dim) const
 Return the coordinate (weight) of the specified node. For straight-edged geometry, this is simply the physical coordinate of the vertex. For all geometries, this can be understood as a weight on the corresponding H^1 basis function used in the reference-to-physical map.
 
KOKKOS_INLINE_FUNCTION int uniformJacobianModulus () const
 Returns an integer indicating the number of distinct cell types vis-a-vis Jacobians.
 
Data< PointScalar, DeviceType > allocateJacobianData (const TensorPoints< PointScalar, DeviceType > &points, const int startCell=0, const int endCell=-1) const
 Allocate a container into which Jacobians of the reference-to-physical mapping can be placed. More...
 
Data< PointScalar, DeviceType > allocateJacobianData (const ScalarView< PointScalar, DeviceType > &points, const int startCell=0, const int endCell=-1) const
 Allocate a container into which Jacobians of the reference-to-physical mapping can be placed. More...
 
Data< PointScalar, DeviceType > allocateJacobianData (const int &numPoints, const int startCell=0, const int endCell=-1) const
 Allocate a container into which Jacobians of the reference-to-physical mapping can be placed (variant for affine geometry). More...
 
Data< PointScalar, DeviceType > getJacobianRefData (const ScalarView< PointScalar, DeviceType > &points) const
 Computes reference-space data for the specified points, to be used in setJacobian(). More...
 
Data< PointScalar, DeviceType > getJacobianRefData (const TensorPoints< PointScalar, DeviceType > &points) const
 Computes reference-space data for the specified points, to be used in setJacobian(). More...
 
void setJacobian (Data< PointScalar, DeviceType > &jacobianData, const TensorPoints< PointScalar, DeviceType > &points, const Data< PointScalar, DeviceType > &refData, const int startCell=0, const int endCell=-1) const
 Compute Jacobian values for the reference-to-physical transformation, and place them in the provided container. More...
 
void setJacobian (Data< PointScalar, DeviceType > &jacobianData, const ScalarView< PointScalar, DeviceType > &points, const Data< PointScalar, DeviceType > &refData, const int startCell=0, const int endCell=-1) const
 Compute Jacobian values for the reference-to-physical transformation, and place them in the provided container. More...
 
void setJacobian (Data< PointScalar, DeviceType > &jacobianData, const int &numPoints, const int startCell=0, const int endCell=-1) const
 Compute Jacobian values for the reference-to-physical transformation, and place them in the provided container (variant for affine geometry). More...
 

Protected Types

using BasisPtr = Teuchos::RCP< Basis< DeviceType, PointScalar, PointScalar > >
 

Protected Attributes

HypercubeNodeOrdering nodeOrdering_
 
CellGeometryType cellGeometryType_
 
SubdivisionStrategy subdivisionStrategy_ = NO_SUBDIVISION
 
bool affine_
 
Data< Orientation, DeviceType > orientations_
 
Kokkos::Array< PointScalar,
spaceDim > 
origin_
 
Kokkos::Array< PointScalar,
spaceDim > 
domainExtents_
 
Kokkos::Array< int, spaceDim > gridCellCounts_
 
TensorPoints< PointScalar,
DeviceType > 
tensorVertices_
 
ScalarView< int, DeviceType > cellToNodes_
 
ScalarView< PointScalar,
DeviceType > 
nodes_
 
unsigned numCells_ = 0
 
unsigned numNodesPerCell_ = 0
 

Detailed Description

template<class PointScalar, int spaceDim, typename DeviceType>
class Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >

CellGeometry provides the nodes for a set of cells; has options that support efficient definition of uniform grids as well as options for arbitrary geometry, including curvilinear.

CellGeometry implements allocate and set methods for cell Jacobians, their inverses, and their determinants. These are computed and stored in a fashion that avoids redundant computation and storage when CellGeometry parameters imply structure, either because the mapping is affine, or because there is axis-aligned tensorial structure that implies a diagonal Jacobian.

CellGeometry also handles orientations for the cells, accessible through the getOrientation() and getOrientations(). These are backed by lazily-initialized storage of orientations for all cells.

Note
Conceptually, it would make sense to use class inheritance and have different member data for each type of geometry supported. We instead glom all the options together into one multi-modal class; this is basically to avoid certain difficulties with vtables under CUDA.

Definition at line 80 of file Intrepid2_CellGeometry.hpp.

Member Enumeration Documentation

template<class PointScalar, int spaceDim, typename DeviceType>
enum Intrepid2::CellGeometry::CellGeometryType

Supported types for CellGeometry

Enumerator
UNIFORM_GRID 

each grid division has the same dimensions

TENSOR_GRID 

grid expressed as a Cartesian product of 1D grids (could be a Shishkin mesh, e.g.)

EXTRUDED_GRID 

lower-dimensional geometry that is orthogonally extruded in higher dimensions

FIRST_ORDER 

geometry expressible in terms of vertices of the cell

HIGHER_ORDER 

geometry expressible in terms of a higher-order basis (must be specified)

Definition at line 84 of file Intrepid2_CellGeometry.hpp.

template<class PointScalar, int spaceDim, typename DeviceType>
enum Intrepid2::CellGeometry::HypercubeNodeOrdering

Distinguish between "classic" (counterclockwise in 2D, counterclockwise bottom followed by counterclockwise top in 3D) Shards ordering and a more natural tensor ordering. The tensor ordering is such that the lowest-order bit in the composite vertex number corresponds to the x dimension, the second-lowest-order bit corresponds to the y dimension, etc. (There is not yet a non-"classic" Shards ordering, but we hope to add one.)

Enumerator
HYPERCUBE_NODE_ORDER_CLASSIC_SHARDS 

classic shards ordering

HYPERCUBE_NODE_ORDER_TENSOR 

a more natural tensor ordering

Definition at line 105 of file Intrepid2_CellGeometry.hpp.

template<class PointScalar, int spaceDim, typename DeviceType>
enum Intrepid2::CellGeometry::SubdivisionStrategy

Options for uniform, tensor grids subdivided into simplices

Enumerator
NO_SUBDIVISION 

no subdivision

TWO_TRIANGLES_RIGHT 

square –> two triangles, with a hypotenuse of slope 1

TWO_TRIANGLES_LEFT 

square –> two triangles, with a hypotenuse of slope -1

FOUR_TRIANGLES 

square –> four triangles, with a new vertex at center

FIVE_TETRAHEDRA 

cube –> five tetrahedra

SIX_TETRAHEDRA 

cube –> six tetrahedra

SIX_PYRAMIDS 

cube –> six pyramids

Definition at line 93 of file Intrepid2_CellGeometry.hpp.

Constructor & Destructor Documentation

template<class PointScalar, int spaceDim, typename DeviceType >
Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::CellGeometry ( const Kokkos::Array< PointScalar, spaceDim > &  origin,
const Kokkos::Array< PointScalar, spaceDim > &  domainExtents,
const Kokkos::Array< int, spaceDim > &  gridCellCounts,
SubdivisionStrategy  subdivisionStrategy = NO_SUBDIVISION,
HypercubeNodeOrdering  nodeOrdering = HYPERCUBE_NODE_ORDER_TENSOR 
)

Uniform grid constructor, with optional subdivision into simplices.

Parameters
[in]origin- location of the "origin" of the mesh, coordinates start here; node locations are based on the sum of origin and domain extents, subdivided according to cell counts in each dimension.
[in]domainExtents- the size of the domain in each coordinate dimension.
[in]gridCellCounts- the number of grid cells in each coordinate dimension.
[in]subdivisionStrategy- whether (and how) to subdivide the (hypercube) grid elements into simplices.
[in]nodeOrdering- applicable for hypercube cell topologies; specifies whether to use the order used by Shards (and lowest-order Intrepid2 bases), or the one used by higher-order Intrepid2 bases.

Definition at line 418 of file Intrepid2_CellGeometryDef.hpp.

References Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::basisForNodes(), Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::HYPERCUBE_NODE_ORDER_CLASSIC_SHARDS, Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::NO_SUBDIVISION, and Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::numCellsPerGridCell().

template<class PointScalar, int spaceDim, typename DeviceType>
Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::CellGeometry ( const shards::CellTopology &  cellTopo,
ScalarView< int, DeviceType >  cellToNodes,
ScalarView< PointScalar, DeviceType >  nodes,
const bool  claimAffine = false,
const HypercubeNodeOrdering  nodeOrdering = HYPERCUBE_NODE_ORDER_CLASSIC_SHARDS 
)

Node-based constructor for straight-edged geometry.

Parameters
[in]cellTopo- the cell topology for all cells.
[in]cellToNodes- (C,LN) container specifying the global index of the local node.
[in]nodes- (GN,D) container specifying the coordinate weight for the global node in the specified dimension; if cellToNodes is not allocated, this must be a (C,N,D) container
[in]claimAffine- whether to assume (without checking) that the mapping from reference space is affine. (If claimAffine is false, we check whether the cell topology is simplicial, and if so, set the affine_ member variable to true.)
[in]nodeOrdering- applicable for hypercube cell topologies; specifies whether to use the order used by Shards (and lowest-order Intrepid2 bases), or the one used by higher-order Intrepid2 bases.

Definition at line 506 of file Intrepid2_CellGeometryDef.hpp.

References Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::basisForNodes(), Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::HYPERCUBE_NODE_ORDER_CLASSIC_SHARDS, and INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE.

template<class PointScalar, int spaceDim, typename DeviceType>
Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::CellGeometry ( Teuchos::RCP< Intrepid2::Basis< DeviceType, PointScalar, PointScalar > >  basisForNodes,
ScalarView< PointScalar, DeviceType >  cellNodes 
)

Node-based constructor for curvilinear geometry.

Parameters
[in]basisForNodes- the basis in terms of which the reference-to-physical transformation is expressed.
[in]cellNodes- (C,F,D) container specifying the coordinate weight of the node; F is the (local) field ordinal for the basis in terms of which the reference-to-physical transformation is expressed.

Definition at line 566 of file Intrepid2_CellGeometryDef.hpp.

References Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::basisForNodes(), Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::FIRST_ORDER, and Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::HIGHER_ORDER.

template<class PointScalar, int spaceDim, typename DeviceType>
KOKKOS_INLINE_FUNCTION Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::CellGeometry ( const CellGeometry< PointScalar, spaceDim, DeviceType > &  cellGeometry)

Member Function Documentation

template<class PointScalar, int spaceDim, typename DeviceType>
TensorData< PointScalar, DeviceType > Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::allocateCellMeasure ( const Data< PointScalar, DeviceType > &  jacobianDet,
const TensorData< PointScalar, DeviceType > &  cubatureWeights 
) const

Allocate a TensorData object appropriate for passing to computeCellMeasure().

Parameters
[in]jacobianDet- a container approriately sized to store the determinants of the Jacobians of the reference-to-physical mapping.
[in]cubatureWeights- a container approriately sized to store the quadrature weights.
Returns
TensorData object sized to accept the result of computeCellMeasure() when called with the provided jacobianDet and cubatureWeights arguments.

Definition at line 603 of file Intrepid2_CellGeometryDef.hpp.

References Intrepid2::TensorData< Scalar, DeviceType >::extent_int(), Intrepid2::Data< DataScalar, DeviceType >::extent_int(), Intrepid2::Data< DataScalar, DeviceType >::getDataExtent(), Intrepid2::TensorData< Scalar, DeviceType >::getTensorComponent(), Intrepid2::Data< DataScalar, DeviceType >::getVariationTypes(), Intrepid2::TensorData< Scalar, DeviceType >::numTensorComponents(), and Intrepid2::TensorData< Scalar, DeviceType >::rank().

template<class PointScalar, int spaceDim, typename DeviceType>
Data< PointScalar, DeviceType > Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::allocateJacobianData ( const TensorPoints< PointScalar, DeviceType > &  points,
const int  startCell = 0,
const int  endCell = -1 
) const

Allocate a container into which Jacobians of the reference-to-physical mapping can be placed.

Parameters
[in]points- the points at which the Jacobian will be evaluated.
[in]startCell- the first cell ordinal for which the Jacobian will be evaluated (used to define worksets smaller than the whole CellGeometry).
[in]endCell- the first cell ordinal for which the Jacobian will not be evaluated (used to define worksets smaller than the whole CellGeometry); use -1 to indicate that all cells from startCell on should be evaluated.
Returns
a Data object appropriately sized to accomodate the specified Jacobian values.

Definition at line 1504 of file Intrepid2_CellGeometryDef.hpp.

References Intrepid2::TensorPoints< PointScalar, DeviceType >::extent_int(), and Intrepid2::TensorPoints< PointScalar, DeviceType >::getTensorComponent().

Referenced by Intrepid2::ProjectedGeometry< spaceDim, PointScalar, DeviceType >::projectOntoHGRADBasis().

template<class PointScalar, int spaceDim, typename DeviceType>
Data< PointScalar, DeviceType > Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::allocateJacobianData ( const ScalarView< PointScalar, DeviceType > &  points,
const int  startCell = 0,
const int  endCell = -1 
) const

Allocate a container into which Jacobians of the reference-to-physical mapping can be placed.

Parameters
[in]points- the points at which the Jacobian will be evaluated.
[in]startCell- the first cell ordinal for which the Jacobian will be evaluated (used to define worksets smaller than the whole CellGeometry).
[in]endCell- the first cell ordinal for which the Jacobian will not be evaluated (used to define worksets smaller than the whole CellGeometry); use -1 to indicate that all cells from startCell on should be evaluated.
Returns
a Data object appropriately sized to accomodate the specified Jacobian values.

Definition at line 1511 of file Intrepid2_CellGeometryDef.hpp.

template<class PointScalar, int spaceDim, typename DeviceType>
Data< PointScalar, DeviceType > Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::allocateJacobianData ( const int &  numPoints,
const int  startCell = 0,
const int  endCell = -1 
) const

Allocate a container into which Jacobians of the reference-to-physical mapping can be placed (variant for affine geometry).

Parameters
[in]numPoints- the number of points at which the Jacobian will be defined.
[in]startCell- the first cell ordinal for which the Jacobian will be evaluated (used to define worksets smaller than the whole CellGeometry).
[in]endCell- the first cell ordinal for which the Jacobian will not be evaluated (used to define worksets smaller than the whole CellGeometry); use -1 to indicate that all cells from startCell on should be evaluated.
Returns
a Data object appropriately sized to accomodate the specified Jacobian values.

Definition at line 1520 of file Intrepid2_CellGeometryDef.hpp.

References INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE.

template<class PointScalar, int spaceDim, typename DeviceType>
Data< PointScalar, DeviceType > Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::allocateJacobianDataPrivate ( const ScalarView< PointScalar, DeviceType > &  pointComponentView,
const int &  pointsPerCell,
const int  startCell,
const int  endCell 
) const

Notionally-private method that provides a common interface for multiple public-facing allocateJacobianData() methods. (Marked as public due to compiler constraints.)

Parameters
[in]pointComponentView- typically either the first component of a TensorPoints container, or a View containing all the points, but may be empty. Only used for scalar-type-size-matched construction of new views (important for views of Fad types with hidden dimensions).
[in]pointsPerCell- the number of points at which the Jacobian will be evaluated in each cell. If points is a valid container, pointsPerCell must match its first dimension.
[in]startCell- the first cell ordinal for which the Jacobian will be evaluated (used to define worksets smaller than the whole CellGeometry).
[in]endCell- the first cell ordinal for which the Jacobian will not be evaluated (used to define worksets smaller than the whole CellGeometry); use -1 to indicate that all cells from startCell on should be evaluated.
Returns
a Data object appropriately sized to accomodate the specified Jacobian values.

Definition at line 185 of file Intrepid2_CellGeometryDef.hpp.

References INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE.

template<class PointScalar , int spaceDim, typename DeviceType >
KOKKOS_INLINE_FUNCTION DataVariationType Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::cellVariationType ( ) const

Indicates the type of geometric variation from one cell to the next. If all cells are known to have the same geometry, returns CONSTANT. Returns CONSTANT for uniform grids with no subdivisions, MODULAR for uniform grids with subdivisions, GENERAL for all others.

Definition at line 814 of file Intrepid2_CellGeometryDef.hpp.

template<class PointScalar, int spaceDim, typename DeviceType>
void Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::computeCellMeasure ( TensorData< PointScalar, DeviceType > &  cellMeasure,
const Data< PointScalar, DeviceType > &  jacobianDet,
const TensorData< PointScalar, DeviceType > &  cubatureWeights 
) const

Compute cell measures that correspond to provided Jacobian determinants and.

Parameters
[out]cellMeasure- a container, usually provided by allocateCellMeasure(), to store the cell measures.
[in]jacobianDet- the determinants of the Jacobians of the reference-to-physical mapping.
[in]cubatureWeights- the quadrature weights.

Definition at line 662 of file Intrepid2_CellGeometryDef.hpp.

References Intrepid2::TensorData< Scalar, DeviceType >::extent_int(), Intrepid2::TensorData< Scalar, DeviceType >::getTensorComponent(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView1(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView2(), Intrepid2::Data< DataScalar, DeviceType >::getVariationTypes(), Intrepid2::TensorData< Scalar, DeviceType >::numTensorComponents(), and Intrepid2::TensorData< Scalar, DeviceType >::rank().

template<class PointScalar, int spaceDim, typename DeviceType>
Data< PointScalar, DeviceType > Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::getJacobianRefData ( const ScalarView< PointScalar, DeviceType > &  points) const

Computes reference-space data for the specified points, to be used in setJacobian().

Parameters
[in]points- the points at which the Jacobian will be evaluated; shape (P,D) or (C,P,D).
Returns
a Data object with any reference-space data required. This may be empty, if no reference-space data is required in setJacobian(). If filled, will have shape (F,P,D) or (C,F,P,D).

Definition at line 832 of file Intrepid2_CellGeometryDef.hpp.

References INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE.

Referenced by Intrepid2::ProjectedGeometry< spaceDim, PointScalar, DeviceType >::projectOntoHGRADBasis().

template<class PointScalar, int spaceDim, typename DeviceType>
Data< PointScalar, DeviceType > Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::getJacobianRefData ( const TensorPoints< PointScalar, DeviceType > &  points) const

Computes reference-space data for the specified points, to be used in setJacobian().

Parameters
[in]points- the points at which the Jacobian will be evaluated, with shape (P,D).
Returns
a Data object with any reference-space data required. This may be empty, if no reference-space data is required in setJacobian(). If filled, will have shape (F,P,D).

Definition at line 916 of file Intrepid2_CellGeometryDef.hpp.

References Intrepid2::TensorPoints< PointScalar, DeviceType >::extent_int(), Intrepid2::TensorPoints< PointScalar, DeviceType >::getTensorComponent(), and INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE.

template<class PointScalar , int spaceDim, typename DeviceType >
KOKKOS_INLINE_FUNCTION int Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::numCellsPerGridCell ( SubdivisionStrategy  subdivisionStrategy) const

Helper method that returns the number of cells into which each grid cell will be subdivided based on the input subdivision strategy.

Parameters
[in]subdivisionStrategy- the subdivision strategy.
Returns
the number of subdivisions specified by the subdivision strategy; -1 for unsupported subdivisionStrategy values.

Definition at line 164 of file Intrepid2_CellGeometryDef.hpp.

Referenced by Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::CellGeometry().

template<class PointScalar , int spaceDim, typename DeviceType>
void Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::orientations ( ScalarView< Orientation, DeviceType >  orientationsView,
const int &  startCell = 0,
const int &  endCell = -1 
)

Fills the provided container with the orientations for the specified cell range. Calls getOrientations() and copies the orientations from the Data container into the ScalarView container.

Parameters
[out]orientationsView- the container that will be filled.
[in]startCell- the first cell ordinal whose orientation will be copied.
[in]endCell- the first cell ordinal whose orientation will not be copied; use -1 to indicate that orientations for all cells from startCell on should be copied.

Definition at line 1475 of file Intrepid2_CellGeometryDef.hpp.

template<class PointScalar, int spaceDim, typename DeviceType>
void Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::setJacobian ( Data< PointScalar, DeviceType > &  jacobianData,
const TensorPoints< PointScalar, DeviceType > &  points,
const Data< PointScalar, DeviceType > &  refData,
const int  startCell = 0,
const int  endCell = -1 
) const

Compute Jacobian values for the reference-to-physical transformation, and place them in the provided container.

Parameters
[out]jacobianData- a container, allocated by allocateJacobianData(), into which the evaluated Jacobians will be placed.
[in]points- the points at which the Jacobian will be evaluated.
[in]refData- the return from getJacobianRefData(); may be an empty container, depending on details of CellGeometry (e.g. if it is affine)
[in]startCell- the first cell ordinal for which the Jacobian will be evaluated (used to define worksets smaller than the whole CellGeometry).
[in]endCell- the first cell ordinal for which the Jacobian will not be evaluated (used to define worksets smaller than the whole CellGeometry); use -1 to indicate that all cells from startCell on should be evaluated.

Definition at line 1529 of file Intrepid2_CellGeometryDef.hpp.

References Intrepid2::TensorPoints< PointScalar, DeviceType >::extent_int().

Referenced by Intrepid2::ProjectedGeometry< spaceDim, PointScalar, DeviceType >::projectOntoHGRADBasis().

template<class PointScalar, int spaceDim, typename DeviceType>
void Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::setJacobian ( Data< PointScalar, DeviceType > &  jacobianData,
const ScalarView< PointScalar, DeviceType > &  points,
const Data< PointScalar, DeviceType > &  refData,
const int  startCell = 0,
const int  endCell = -1 
) const

Compute Jacobian values for the reference-to-physical transformation, and place them in the provided container.

Parameters
[out]jacobianData- a container, allocated by allocateJacobianData(), into which the evaluated Jacobians will be placed.
[in]points- the points at which the Jacobian will be evaluated.
[in]refData- the return from getJacobianRefData(); may be an empty container, depending on details of CellGeometry (e.g. if it is affine)
[in]startCell- the first cell ordinal for which the Jacobian will be evaluated (used to define worksets smaller than the whole CellGeometry).
[in]endCell- the first cell ordinal for which the Jacobian will not be evaluated (used to define worksets smaller than the whole CellGeometry); use -1 to indicate that all cells from startCell on should be evaluated.

Definition at line 1537 of file Intrepid2_CellGeometryDef.hpp.

template<class PointScalar, int spaceDim, typename DeviceType>
void Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::setJacobian ( Data< PointScalar, DeviceType > &  jacobianData,
const int &  numPoints,
const int  startCell = 0,
const int  endCell = -1 
) const

Compute Jacobian values for the reference-to-physical transformation, and place them in the provided container (variant for affine geometry).

Parameters
[out]jacobianData- a container, allocated by allocateJacobianData(), into which the evaluated Jacobians will be placed.
[in]numPoints- the number of points at which the Jacobian will be defined.
[in]startCell- the first cell ordinal for which the Jacobian will be evaluated (used to define worksets smaller than the whole CellGeometry).
[in]endCell- the first cell ordinal for which the Jacobian will not be evaluated (used to define worksets smaller than the whole CellGeometry); use -1 to indicate that all cells from startCell on should be evaluated.

Definition at line 1547 of file Intrepid2_CellGeometryDef.hpp.

References INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE.

template<class PointScalar, int spaceDim, typename DeviceType>
void Intrepid2::CellGeometry< PointScalar, spaceDim, DeviceType >::setJacobianDataPrivate ( Data< PointScalar, DeviceType > &  jacobianData,
const int &  pointsPerCell,
const Data< PointScalar, DeviceType > &  refData,
const int  startCell,
const int  endCell 
) const

Notionally-private method that provides a common interface for multiple public-facing setJacobianData() methods. (Marked as public due to compiler constraints.)

Parameters
[out]jacobianData- a container, allocated by allocateJacobianData(), into which the evaluated Jacobians will be placed.
[in]pointsPerCell- the number of points at which the Jacobian will be evaluated in each cell. If points is a valid container, pointsPerCell must match its first dimension.
[in]refData- the return from getJacobianRefData(); may be an empty container, depending on details of CellGeometry (e.g. if it is affine)
[in]startCell- the first cell ordinal for which the Jacobian will be evaluated (used to define worksets smaller than the whole CellGeometry).
[in]endCell- the first cell ordinal for which the Jacobian will not be evaluated (used to define worksets smaller than the whole CellGeometry); use -1 to indicate that all cells from startCell on should be evaluated.

Definition at line 277 of file Intrepid2_CellGeometryDef.hpp.

References Intrepid2::Data< DataScalar, DeviceType >::extent_int(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView1(), Intrepid2::Data< DataScalar, DeviceType >::getUnderlyingView3(), INTREPID2_TEST_FOR_EXCEPTION_DEVICE_SAFE, Intrepid2::Data< DataScalar, DeviceType >::isValid(), Intrepid2::Data< DataScalar, DeviceType >::rank(), and Intrepid2::CellTools< DeviceType >::setJacobian().


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