19#ifndef __INTREPID2_CELLDATA_HPP__
20#define __INTREPID2_CELLDATA_HPP__
22#include "Intrepid2_ConfigDefs.hpp"
24#include "Shards_CellTopology.hpp"
76template<
typename DeviceType>
79 using ConstViewType = Kokkos::DynRankView<const double,DeviceType>;
115 get(
const ordinal_type subcellDim,
116 const unsigned parentCellKey );
134 template <
typename HostViewType>
136 set( HostViewType subcellParam,
137 const ordinal_type subcellDim,
138 const shards::CellTopology parentCell );
141 using ViewType = Kokkos::DynRankView<double,DeviceType>;
143 static ViewType triEdgesParam, quadEdgesParam;
144 static ViewType shellTriEdgesParam, shellQuadEdgesParam;
145 static ViewType tetEdgesParam, hexEdgesParam, pyrEdgesParam, wedgeEdgesParam;
146 static ViewType shellTriFacesParam, shellQuadFacesParam;
147 static ViewType tetFacesParam, hexFacesParam, pyrFacesParam, wedgeFacesParam;
171template<
typename DeviceType>
174 using ConstViewType = Kokkos::DynRankView<const double,DeviceType>;
200 get(
const unsigned cellTopoKey);
208 using ViewType = Kokkos::DynRankView<double,DeviceType>;
209 static ViewType lineNodes, line3Nodes;
210 static ViewType triangleNodes, triangle4Nodes, triangle6Nodes;
211 static ViewType quadrilateralNodes, quadrilateral8Nodes, quadrilateral9Nodes;
212 static ViewType tetrahedronNodes, tetrahedron8Nodes, tetrahedron10Nodes, tetrahedron11Nodes;
213 static ViewType hexahedronNodes, hexahedron20Nodes, hexahedron27Nodes;
214 static ViewType pyramidNodes, pyramid13Nodes, pyramid14Nodes;
215 static ViewType wedgeNodes, wedge15Nodes, wedge18Nodes;
222 double line[2][3], line_3[3][3];
223 double triangle[3][3], triangle_4[4][3], triangle_6[6][3];
224 double quadrilateral[4][3], quadrilateral_8[8][3], quadrilateral_9[9][3];
225 double tetrahedron[4][3], tetrahedron_8[8][3], tetrahedron_10[10][3], tetrahedron_11[10][3];
226 double hexahedron[8][3], hexahedron_20[20][3], hexahedron_27[27][3];
227 double pyramid[5][3], pyramid_13[13][3], pyramid_14[14][3];
228 double wedge[6][3], wedge_15[15][3], wedge_18[18][3];
253template<
typename DeviceType>
256 using ConstViewType = Kokkos::DynRankView<const double,DeviceType>;
282 get(
const unsigned cellTopoKey);
290 using ViewType = Kokkos::DynRankView<double,DeviceType>;
293 static ViewType quadrilateralCenter;
305 double quadrilateral[3];
306 double tetrahedron[3];
307 double hexahedron[3];
332template<
unsigned CellTopologyKey>
340 template<
typename Po
intViewType,
typename ScalarType>
341 KOKKOS_INLINE_FUNCTION
343 check(
const PointViewType &point,
const ScalarType threshold);
351 template<
typename Po
intViewType,
typename ScalarType>
352 KOKKOS_INLINE_FUNCTION
354 check(
const PointViewType &point,
const ScalarType threshold);
363 template<
typename Po
intViewType,
typename ScalarType>
364 KOKKOS_INLINE_FUNCTION
366 check(
const PointViewType &point,
const ScalarType threshold);
374 template<
typename Po
intViewType,
typename ScalarType>
375 KOKKOS_INLINE_FUNCTION
377 check(
const PointViewType &point,
const ScalarType threshold);
385 template<
typename Po
intViewType,
typename ScalarType>
386 KOKKOS_INLINE_FUNCTION
388 check(
const PointViewType &point,
const ScalarType threshold);
396 template<
typename Po
intViewType,
typename ScalarType>
397 KOKKOS_INLINE_FUNCTION
399 check(
const PointViewType &point,
const ScalarType threshold);
407 template<
typename Po
intViewType,
typename ScalarType>
408 KOKKOS_INLINE_FUNCTION
410 check(
const PointViewType &point,
const ScalarType threshold);
413 const CellTopologyData* getCellTopologyData(
const unsigned& cellTopologyKey);
Definition file for the classes: Intrepid2::RefSubcellParametrization, Intrepid2::RefCellNodes,...
Header file for small functions used in Intrepid2.
Contains definitions of custom data types in Intrepid2.
Header function for Intrepid2::Util class and other utility functions.
This class defines the coordinates of the barycenter of the supported reference cells....
~RefCellCenter()=default
Destructor.
Kokkos::DynRankView< double, DeviceType > ViewType
static views containing the center coordinates allocated on DeviceType::memory_space
static const ReferenceCenterDataStatic refCenterDataStatic_
static struct containing the nodes coordinates on host
static void set()
Set center coordinates of reference cell for supported topologies.
RefCellCenter()=default
Default constructor.
static ConstViewType get(const unsigned cellTopoKey)
Retrieves the Cartesian coordinates of a reference cell barycenter.
static bool isReferenceCellCenterDataSet_
whether the center coordinates have been already set using the method set()
This class defines the coordinates of the nodes of reference cells according for supported cell topol...
~RefCellNodes()=default
Destructor.
static ConstViewType get(const unsigned cellTopoKey)
Retrieves the Cartesian coordinates of reference cell nodes.
static void set()
Set reference nodes coordinates for supported topologies.
static bool isReferenceNodeDataSet_
whether the nodes coordinates have been already set using the method set()
static const ReferenceNodeDataStatic refNodeDataStatic_
static struct containing the nodes coordinates on host
RefCellNodes()=default
Default constructor.
Kokkos::DynRankView< double, DeviceType > ViewType
static views containing the node coordinates allocated on DeviceType::memory_space
This class defines the parametrizations of edges and faces of supported reference cells....
static void set()
Computes and stores static views containing the parametrizations maps of edges and faces of all refer...
static bool isSubcellParametrizationSet_
whether the parametrizations have been already computed using the method set()
~RefSubcellParametrization()=default
Destructor.
Kokkos::DynRankView< double, DeviceType > ViewType
static views containing the parametrization maps, allocated on DeviceType::memory_space
static bool isSupported(const unsigned cellTopoKey)
Checks if a cell topology has a reference parametrization.
RefSubcellParametrization()=default
Default constructor.
static ConstViewType get(const ordinal_type subcellDim, const unsigned parentCellKey)
Returns a Kokkos view with the coefficients of the parametrization maps for the edges or faces of a r...
This class implements a check function that determines whether a given point is inside or outside the...
Reference node containers for each supported topology.