Intrepid2
Intrepid2_CellGeometry.hpp
Go to the documentation of this file.
1// @HEADER
2// ************************************************************************
3//
4// Intrepid2 Package
5// Copyright (2007) Sandia Corporation
6//
7// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8// license for use of this work by or on behalf of the U.S. Government.
9//
10// Redistribution and use in source and binary forms, with or without
11// modification, are permitted provided that the following conditions are
12// met:
13//
14// 1. Redistributions of source code must retain the above copyright
15// notice, this list of conditions and the following disclaimer.
16//
17// 2. Redistributions in binary form must reproduce the above copyright
18// notice, this list of conditions and the following disclaimer in the
19// documentation and/or other materials provided with the distribution.
20//
21// 3. Neither the name of the Corporation nor the names of the
22// contributors may be used to endorse or promote products derived from
23// this software without specific prior written permission.
24//
25// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36//
37// Questions? Contact Kyungjoo Kim (kyukim@sandia.gov),
38// Mauro Perego (mperego@sandia.gov), or
39// Nate Roberts (nvrober@sandia.gov),
40//
41// ************************************************************************
42// @HEADER
43
50#ifndef Intrepid2_CellGeometry_h
51#define Intrepid2_CellGeometry_h
52
53#include "Intrepid2_Basis.hpp"
55#include "Intrepid2_Data.hpp"
62#include "Intrepid2_Utils.hpp"
64
65#include "Intrepid2_ScalarView.hpp"
66
67namespace Intrepid2
68{
79 template<class PointScalar, int spaceDim, typename DeviceType>
81 {
82 public:
90 };
91
94 {
101 };
102
105 {
108 };
109
114 KOKKOS_INLINE_FUNCTION
115 int numCellsPerGridCell(SubdivisionStrategy subdivisionStrategy) const;
116
117 public:
125 Data<PointScalar,DeviceType> allocateJacobianDataPrivate(const TensorPoints<PointScalar,DeviceType> &points, const int &pointsPerCell, const int startCell, const int endCell) const;
126
135 void setJacobianDataPrivate(Data<PointScalar,DeviceType> &jacobianData, const TensorPoints<PointScalar,DeviceType> &points, const int &pointsPerCell,
136 const Data<PointScalar,DeviceType> &refData, const int startCell, const int endCell) const;
137 protected:
138 HypercubeNodeOrdering nodeOrdering_;
139 CellGeometryType cellGeometryType_;
140 SubdivisionStrategy subdivisionStrategy_ = NO_SUBDIVISION;
141 bool affine_; // if true, each cell has constant Jacobian across the cell
142 Data<Orientation, DeviceType> orientations_; // for grid types, this could have either a single entry or one matching numCellsPerGridCell(). For other types, it has as many entries as there are cells.
143
144 // uniform grid data -- used for UNIFORM_GRID type
145 Kokkos::Array<PointScalar,spaceDim> origin_; // point specifying a corner of the mesh
146 Kokkos::Array<PointScalar,spaceDim> domainExtents_; // how far the domain extends in each dimension
147 Kokkos::Array<int,spaceDim> gridCellCounts_; // how many grid cells wide the mesh is in each dimension
148
149 // tensor grid data -- only used for TENSOR_GRID type
151
152 // arbitrary cell node data, used for both higher-order and first-order
153 // (here, nodes are understood as geometry degrees of freedom)
154 ScalarView<int,DeviceType> cellToNodes_; // (C,N) -- N is the number of nodes per cell; values are global node ordinals
155 ScalarView<PointScalar,DeviceType> nodes_; // (GN,D) or (C,N,D) -- GN is the number of global nodes; (C,N,D) used only if cellToNodes_ is empty.
156 using BasisPtr = Teuchos::RCP< Basis<DeviceType,PointScalar,PointScalar> >;
157
158 unsigned numCells_ = 0;
159 unsigned numNodesPerCell_ = 0;
160 public:
168 CellGeometry(const Kokkos::Array<PointScalar,spaceDim> &origin,
169 const Kokkos::Array<PointScalar,spaceDim> &domainExtents,
170 const Kokkos::Array<int,spaceDim> &gridCellCounts,
171 SubdivisionStrategy subdivisionStrategy = NO_SUBDIVISION,
173
181 CellGeometry(const shards::CellTopology &cellTopo,
182 ScalarView<int,DeviceType> cellToNodes,
183 ScalarView<PointScalar,DeviceType> nodes,
184 const bool claimAffine = false,
186
192 ScalarView<PointScalar,DeviceType> cellNodes);
193
197 KOKKOS_INLINE_FUNCTION CellGeometry(const CellGeometry &cellGeometry);
198
201 KOKKOS_INLINE_FUNCTION ~CellGeometry();
202
204 KOKKOS_INLINE_FUNCTION
205 bool affine() const;
206
213
219 void computeCellMeasure( TensorData<PointScalar,DeviceType> &cellMeasure, const Data<PointScalar,DeviceType> & jacobianDet, const TensorData<PointScalar,DeviceType> & cubatureWeights ) const;
220
222 BasisPtr basisForNodes() const;
223
225 const shards::CellTopology & cellTopology() const;
226
229 KOKKOS_INLINE_FUNCTION
231
233 KOKKOS_INLINE_FUNCTION
234 int hypercubeComponentNodeNumber(int hypercubeNodeNumber, int d) const;
235
238
240 KOKKOS_INLINE_FUNCTION
241 size_t extent(const int& r) const;
242
244 template <typename iType>
245 KOKKOS_INLINE_FUNCTION
246 typename std::enable_if<std::is_integral<iType>::value, int>::type
247 extent_int(const iType& r) const;
248
250 KOKKOS_INLINE_FUNCTION
252
254 KOKKOS_INLINE_FUNCTION
255 int numCells() const;
256
258 KOKKOS_INLINE_FUNCTION
259 int numCellsInDimension(const int &dim) const;
260
262 KOKKOS_INLINE_FUNCTION
263 int numNodesPerCell() const;
264
266 KOKKOS_INLINE_FUNCTION
267 Orientation getOrientation(int &cellNumber) const;
268
271
273 KOKKOS_INLINE_FUNCTION
274 PointScalar gridCellCoordinate(const int &gridCellOrdinal, const int &localNodeNumber, const int &dim) const;
275
277 KOKKOS_INLINE_FUNCTION
278 unsigned rank() const;
279
281 KOKKOS_INLINE_FUNCTION
282 int gridCellNodeForSubdivisionNode(const int &gridCellOrdinal, const int &subdivisionOrdinal,
283 const int &subdivisionNodeNumber) const;
284
286 KOKKOS_INLINE_FUNCTION
287 PointScalar subdivisionCoordinate(const int &gridCellOrdinal, const int &subdivisionOrdinal,
288 const int &subdivisionNodeNumber, const int &d) const;
289
291 KOKKOS_INLINE_FUNCTION
292 PointScalar
293 operator()(const int& cell, const int& node, const int& dim) const;
294
296 KOKKOS_INLINE_FUNCTION
297 int uniformJacobianModulus() const;
298
305 Data<PointScalar,DeviceType> allocateJacobianData(const TensorPoints<PointScalar,DeviceType> &points, const int startCell=0, const int endCell=-1) const;
306
313 Data<PointScalar,DeviceType> allocateJacobianData(const ScalarView<PointScalar,DeviceType> &points, const int startCell=0, const int endCell=-1) const;
314
321 Data<PointScalar,DeviceType> allocateJacobianData(const int &numPoints, const int startCell=0, const int endCell=-1) const;
322
328
337 const int startCell=0, const int endCell=-1) const;
338
346 void setJacobian(Data<PointScalar,DeviceType> &jacobianData, const ScalarView<PointScalar,DeviceType> &points, const Data<PointScalar,DeviceType> &refData,
347 const int startCell=0, const int endCell=-1) const;
348
355 void setJacobian(Data<PointScalar,DeviceType> &jacobianData, const int &numPoints, const int startCell=0, const int endCell=-1) const;
356 };
357} // namespace Intrepid2
358
359#include <Intrepid2_CellGeometryDef.hpp>
360
361#endif /* Intrepid2_CellGeometry_h */
Header file for the abstract base class Intrepid2::Basis.
Teuchos::RCP< Basis< DeviceType, OutputType, PointType > > BasisPtr
Basis Pointer.
Header file for the Intrepid2::CellTools class.
Defines the Data class, a wrapper around a Kokkos::View that allows data that is constant or repeatin...
Stateless class that acts as a factory for a family of nodal bases (hypercube topologies only at this...
Header file for the Intrepid2::OrientationTools and Intrepid2::Impl::OrientationTools classes.
Header file for the Intrepid2::Orientation class.
View-like interface to tensor data; tensor components are stored separately and multiplied together a...
View-like interface to tensor points; point components are stored separately; the appropriate coordin...
Structure-preserving representation of transformed vector data; reference space values and transforma...
Header function for Intrepid2::Util class and other utility functions.
Reference-space field values for a basis, designed to support typical vector-valued bases.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
CellGeometry provides the nodes for a set of cells; has options that support efficient definition of ...
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.
Data< PointScalar, DeviceType > allocateJacobianDataPrivate(const TensorPoints< PointScalar, DeviceType > &points, const int &pointsPerCell, const int startCell, const int endCell) const
Notionally-private method that provides a common interface for multiple public-facing allocateJacobia...
BasisPtr basisForNodes() const
H^1 Basis used in the reference-to-physical transformation. Linear for straight-edged geometry; highe...
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 ...
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 ...
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 ...
TensorData< PointScalar, DeviceType > allocateCellMeasure(const Data< PointScalar, DeviceType > &jacobianDet, const TensorData< PointScalar, DeviceType > &cubatureWeights) const
Allocate a TensorData object appropriate for passing to computeCellMeasure().
KOKKOS_INLINE_FUNCTION int numCells() const
Returns the number of cells.
KOKKOS_INLINE_FUNCTION int numNodesPerCell() const
Returns the number of nodes per cell; may be more than the number of vertices in the corresponding Ce...
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...
KOKKOS_INLINE_FUNCTION Orientation getOrientation(int &cellNumber) const
Returns the orientation for the specified cell. Requires that initializeOrientations() has been calle...
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
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.
void initializeOrientations()
Initialize the internal orientations_ member with the orientations of each member cell....
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 dimen...
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
@ TWO_TRIANGLES_RIGHT
square --> two triangles, with a hypotenuse of slope 1
@ SIX_TETRAHEDRA
cube --> six tetrahedra
@ FIVE_TETRAHEDRA
cube --> five tetrahedra
@ TWO_TRIANGLES_LEFT
square --> two triangles, with a hypotenuse of slope -1
@ FOUR_TRIANGLES
square --> four triangles, with a new vertex at center
KOKKOS_INLINE_FUNCTION unsigned rank() const
Returns the logical rank of this container. This is always 3.
@ UNIFORM_GRID
each grid division has the same 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)
@ EXTRUDED_GRID
lower-dimensional geometry that is orthogonally extruded in higher dimensions
@ TENSOR_GRID
grid expressed as a Cartesian product of 1D grids (could be a Shishkin mesh, e.g.)
KOKKOS_INLINE_FUNCTION bool affine() const
Returns true if Jacobian is constant within each cell.
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.
Data< PointScalar, DeviceType > getJacobianRefData(const TensorPoints< PointScalar, DeviceType > &points) const
Computes reference-space data for the specified points, to be used in setJacobian().
@ HYPERCUBE_NODE_ORDER_CLASSIC_SHARDS
classic shards ordering
@ HYPERCUBE_NODE_ORDER_TENSOR
a more natural tensor ordering
KOKKOS_INLINE_FUNCTION ~CellGeometry()
Destructor.
KOKKOS_INLINE_FUNCTION HypercubeNodeOrdering nodeOrderingForHypercubes() const
Returns the node ordering used for hypercubes.
Data< Orientation, DeviceType > getOrientations()
Returns the orientations for all cells. Calls initializeOrientations() if it has not previously been ...
void setJacobianDataPrivate(Data< PointScalar, DeviceType > &jacobianData, const TensorPoints< PointScalar, DeviceType > &points, 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...
KOKKOS_INLINE_FUNCTION int uniformJacobianModulus() const
Returns an integer indicating the number of distinct cell types vis-a-vis Jacobians.
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
const shards::CellTopology & cellTopology() const
The shards CellTopology for each cell within the CellGeometry object. Note that this is always a lowe...
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,...
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 CellGe...
Wrapper around a Kokkos::View that allows data that is constant or repeating in various logical dimen...
Orientation encoding and decoding.
View-like interface to tensor points; point components are stored separately; the appropriate coordin...