49 #ifndef __INTREPID2_CELLTOOLS_HPP__
50 #define __INTREPID2_CELLTOOLS_HPP__
52 #include "Intrepid2_ConfigDefs.hpp"
54 #include "Shards_CellTopology.hpp"
55 #include "Shards_BasicTopologies.hpp"
57 #include "Teuchos_RCP.hpp"
109 template<
typename DeviceType>
111 using ExecSpaceType =
typename DeviceType::execution_space;
112 using MemSpaceType =
typename DeviceType::memory_space;
130 template<
typename outputValueType,
131 typename pointValueType>
132 static Teuchos::RCP<Basis<DeviceType,outputValueType,pointValueType> >
134 Teuchos::RCP<Basis<DeviceType,outputValueType,pointValueType> > r_val;
136 switch (cellTopo.getKey()) {
157 case shards::Beam<2>::key:
158 case shards::Beam<3>::key:
159 case shards::ShellLine<2>::key:
160 case shards::ShellLine<3>::key:
161 case shards::ShellTriangle<3>::key:
162 case shards::ShellTriangle<6>::key:
163 case shards::ShellQuadrilateral<4>::key:
164 case shards::ShellQuadrilateral<8>::key:
165 case shards::ShellQuadrilateral<9>::key:
167 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
168 ">>> ERROR (Intrepid2::CellTools::createHGradBasis): Cell topology not supported.");
225 template<
typename JacobianViewType,
226 typename PointViewType,
227 typename WorksetType,
228 typename HGradBasisType>
231 const PointViewType points,
232 const WorksetType worksetCell,
233 const Teuchos::RCP<HGradBasisType> basis,
234 const int startCell=0,
const int endCell=-1);
270 template<
typename JacobianViewType,
271 typename BasisGradientsType,
272 typename WorksetType>
275 const WorksetType worksetCell,
276 const BasisGradientsType gradients,
277 const int startCell=0,
const int endCell=-1);
313 template<
typename JacobianViewType,
314 typename PointViewType,
315 typename WorksetCellViewType>
318 const PointViewType points,
319 const WorksetCellViewType worksetCell,
320 const shards::CellTopology cellTopo ) {
321 using nonConstPointValueType =
typename PointViewType::non_const_value_type;
322 auto basis = createHGradBasis<nonConstPointValueType,nonConstPointValueType>(cellTopo);
339 template<
typename JacobianInvViewType,
340 typename JacobianViewType>
343 const JacobianViewType jacobian );
355 template<
typename JacobianDetViewType,
356 typename JacobianViewType>
359 const JacobianViewType jacobian );
366 template<
class Po
intScalar>
374 template<
class Po
intScalar>
382 template<
class Po
intScalar>
391 template<
class Po
intScalar>
400 template<
class Po
intScalar>
410 template<
class Po
intScalar>
431 template<
typename cellCenterValueType,
class ...cellCenterProperties>
434 const shards::CellTopology cell );
444 template<
typename cellVertexValueType,
class ...cellVertexProperties>
446 getReferenceVertex( Kokkos::DynRankView<cellVertexValueType,cellVertexProperties...> cellVertex,
447 const shards::CellTopology cell,
448 const ordinal_type vertexOrd );
465 template<
typename subcellVertexValueType,
class ...subcellVertexProperties>
468 const ordinal_type subcellDim,
469 const ordinal_type subcellOrd,
470 const shards::CellTopology parentCell );
489 template<
typename cellNodeValueType,
class ...cellNodeProperties>
491 getReferenceNode( Kokkos::DynRankView<cellNodeValueType,cellNodeProperties...> cellNode,
492 const shards::CellTopology cell,
493 const ordinal_type nodeOrd );
510 template<
typename SubcellNodeViewType>
513 const ordinal_type subcellDim,
514 const ordinal_type subcellOrd,
515 const shards::CellTopology parentCell );
542 template<
typename RefEdgeTangentViewType>
545 const ordinal_type edgeOrd,
546 const shards::CellTopology parentCell );
584 template<
typename RefFaceTanViewType>
587 RefFaceTanViewType refFaceTanV,
588 const ordinal_type faceOrd,
589 const shards::CellTopology parentCell );
653 template<
typename RefS
ideNormalViewType>
656 const ordinal_type sideOrd,
657 const shards::CellTopology parentCell );
697 template<
typename RefFaceNormalViewType>
700 const ordinal_type faceOrd,
701 const shards::CellTopology parentCell );
732 template<
typename edgeTangentValueType,
class ...edgeTangentProperties,
733 typename worksetJacobianValueType,
class ...worksetJacobianProperties>
736 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
737 const ordinal_type worksetEdgeOrd,
738 const shards::CellTopology parentCell );
752 template<
typename edgeTangentValueType,
class ...edgeTangentProperties,
753 typename worksetJacobianValueType,
class ...worksetJacobianProperties,
754 typename edgeOrdValueType,
class ...edgeOrdProperties>
757 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
758 const Kokkos::DynRankView<edgeOrdValueType,edgeOrdProperties...> worksetEdgeOrds,
759 const shards::CellTopology parentCell );
800 template<
typename faceTanValueType,
class ...faceTanProperties,
801 typename worksetJacobianValueType,
class ...worksetJacobianProperties>
804 Kokkos::DynRankView<faceTanValueType,faceTanProperties...> faceTanV,
805 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
806 const ordinal_type worksetFaceOrd,
807 const shards::CellTopology parentCell );
822 template<
typename faceTanValueType,
class ...faceTanProperties,
823 typename worksetJacobianValueType,
class ...worksetJacobianProperties,
824 typename faceOrdValueType,
class ...faceOrdProperties>
827 Kokkos::DynRankView<faceTanValueType,faceTanProperties...> faceTanV,
828 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
829 const Kokkos::DynRankView<faceOrdValueType,faceOrdProperties...> worksetFaceOrds,
830 const shards::CellTopology parentCell );
894 template<
typename sideNormalValueType,
class ...sideNormalProperties,
895 typename worksetJacobianValueType,
class ...worksetJacobianProperties>
898 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
899 const ordinal_type worksetSideOrd,
900 const shards::CellTopology parentCell );
913 template<
typename sideNormalValueType,
class ...sideNormalProperties,
914 typename worksetJacobianValueType,
class ...worksetJacobianProperties,
915 typename edgeOrdValueType,
class ...edgeOrdProperties>
918 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
919 const Kokkos::DynRankView<edgeOrdValueType,edgeOrdProperties...> worksetSideOrds,
920 const shards::CellTopology parentCell );
960 template<
typename faceNormalValueType,
class ...faceNormalProperties,
961 typename worksetJacobianValueType,
class ...worksetJacobianProperties>
964 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
965 const ordinal_type worksetFaceOrd,
966 const shards::CellTopology parentCell );
980 template<
typename faceNormalValueType,
class ...faceNormalProperties,
981 typename worksetJacobianValueType,
class ...worksetJacobianProperties,
982 typename faceOrdValueType,
class ...faceOrdProperties>
985 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
986 const Kokkos::DynRankView<faceOrdValueType,faceOrdProperties...> worksetFaceOrds,
987 const shards::CellTopology parentCell );
1031 template<
typename PhysPointValueType,
1032 typename RefPointValueType,
1033 typename WorksetType,
1034 typename HGradBasisPtrType>
1037 const RefPointValueType refPoints,
1038 const WorksetType worksetCell,
1039 const HGradBasisPtrType basis );
1081 template<
typename PhysPointViewType,
1082 typename RefPointViewType,
1083 typename WorksetCellViewType>
1086 const RefPointViewType refPoints,
1087 const WorksetCellViewType worksetCell,
1088 const shards::CellTopology cellTopo ) {
1089 using nonConstRefPointValueType =
typename RefPointViewType::non_const_value_type;
1090 auto basis = createHGradBasis<nonConstRefPointValueType,nonConstRefPointValueType>(cellTopo);
1148 template<
typename refSubcellViewType,
1149 typename paramPointViewType>
1152 const paramPointViewType paramPoints,
1153 const ordinal_type subcellDim,
1154 const ordinal_type subcellOrd,
1155 const shards::CellTopology parentCell );
1164 template<
typename refSubcellViewType,
typename paramPo
intViewType>
1167 const paramPointViewType paramPoints,
1168 const typename RefSubcellParametrization<DeviceType>::ConstViewType subcellParametrization,
1169 const ordinal_type subcellOrd);
1176 template<
typename refSubcellViewType,
typename paramPo
intViewType,
typename ordViewType>
1179 const paramPointViewType paramPoints,
1180 const typename RefSubcellParametrization<DeviceType>::ConstViewType subcellParametrization,
1181 const ordViewType subcellOrd);
1234 template<
typename refPointValueType,
class ...refPointProperties,
1235 typename physPointValueType,
class ...physPointProperties,
1236 typename worksetCellValueType,
class ...worksetCellProperties>
1238 mapToReferenceFrame( Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
1239 const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1240 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1241 const shards::CellTopology cellTopo );
1269 template<
typename refPointValueType,
class ...refPointProperties,
1270 typename initGuessValueType,
class ...initGuessProperties,
1271 typename physPointValueType,
class ...physPointProperties,
1272 typename worksetCellValueType,
class ...worksetCellProperties,
1273 typename HGradBasisPtrType>
1276 const Kokkos::DynRankView<initGuessValueType,initGuessProperties...> initGuess,
1277 const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1278 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1279 const HGradBasisPtrType basis );
1312 template<
typename refPointValueType,
class ...refPointProperties,
1313 typename initGuessValueType,
class ...initGuessProperties,
1314 typename physPointValueType,
class ...physPointProperties,
1315 typename worksetCellValueType,
class ...worksetCellProperties>
1318 const Kokkos::DynRankView<initGuessValueType,initGuessProperties...> initGuess,
1319 const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1320 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1321 const shards::CellTopology cellTopo ) {
1322 auto basis = createHGradBasis<refPointValueType,refPointValueType>(cellTopo);
1410 template<
typename subcvCoordValueType,
class ...subcvCoordProperties,
1411 typename cellCoordValueType,
class ...cellCoordProperties>
1413 getSubcvCoords( Kokkos::DynRankView<subcvCoordValueType,subcvCoordProperties...> subcvCoords,
1414 const Kokkos::DynRankView<cellCoordValueType,cellCoordProperties...> cellCoords,
1415 const shards::CellTopology primaryCell );
1432 template<
typename pointValueType,
class ...pointProperties>
1435 const shards::CellTopology cellTopo,
1436 const double thres = threshold() );
1471 template<
typename inCellValueType,
class ...inCellProperties,
1472 typename pointValueType,
class ...pointProperties>
1474 const Kokkos::DynRankView<pointValueType,pointProperties...> points,
1475 const shards::CellTopology cellTopo,
1476 const double thres = threshold() );
1499 template<
typename inCellValueType,
class ...inCellProperties,
1500 typename pointValueType,
class ...pointProperties,
1501 typename cellWorksetValueType,
class ...cellWorksetProperties>
1503 const Kokkos::DynRankView<pointValueType,pointProperties...> points,
1504 const Kokkos::DynRankView<cellWorksetValueType,cellWorksetProperties...> cellWorkset,
1505 const shards::CellTopology cellTopo,
1506 const double thres = threshold() );
1551 template<
typename jacobianViewType,
1552 typename PointViewType,
1553 typename worksetCellViewType>
1555 CellTools_setJacobianArgs(
const jacobianViewType jacobian,
1556 const PointViewType points,
1557 const worksetCellViewType worksetCell,
1558 const shards::CellTopology cellTopo );
1564 template<
typename jacobianInvViewType,
1565 typename jacobianViewType>
1567 CellTools_setJacobianInvArgs(
const jacobianInvViewType jacobianInv,
1568 const jacobianViewType jacobian );
1575 template<
typename jacobianDetViewType,
1576 typename jacobianViewType>
1578 CellTools_setJacobianDetArgs(
const jacobianDetViewType jacobianDet,
1579 const jacobianViewType jacobian );
1588 template<
typename physPointViewType,
1589 typename refPointViewType,
1590 typename worksetCellViewType>
1592 CellTools_mapToPhysicalFrameArgs(
const physPointViewType physPoints,
1593 const refPointViewType refPoints,
1594 const worksetCellViewType worksetCell,
1595 const shards::CellTopology cellTopo );
1604 template<
typename refPointViewType,
1605 typename physPointViewType,
1606 typename worksetCellViewType>
1608 CellTools_mapToReferenceFrameArgs(
const refPointViewType refPoints,
1609 const physPointViewType physPoints,
1610 const worksetCellViewType worksetCell,
1611 const shards::CellTopology cellTopo );
1622 template<
typename refPointViewType,
1623 typename initGuessViewType,
1624 typename physPointViewType,
1625 typename worksetCellViewType>
1627 CellTools_mapToReferenceFrameInitGuess(
const refPointViewType refPoints,
1628 const initGuessViewType initGuess,
1629 const physPointViewType physPoints,
1630 const worksetCellViewType worksetCell,
1631 const shards::CellTopology cellTopo );
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Quadrilateral cell...
Header file for the Intrepid2::Basis_HGRAD_TRI_C1_FEM class.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Wedge cell.
Implementation of an H(grad)-compatible FEM basis of degree 2 on a Pyramid cell.
Header file for the Intrepid2::Basis_HGRAD_PYR_I2_FEM class.
Header file for the Intrepid2::Basis_HGRAD_QUAD_C1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_QUAD_C2_FEM class.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Triangle cell...
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Pyramid cell...
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Hexahedron cell...
Header file for the Intrepid2::Basis_HGRAD_HEX_C2_FEM class.
Defines the Data class, a wrapper around a Kokkos::View that allows data that is constant or repeatin...
Header file for the Intrepid2::Basis_HGRAD_PYR_C1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_LINE_C2_FEM class.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Triangle cell...
Header file for the Intrepid2::Basis_HGRAD_TRI_C2_FEM class.
Header file for the Intrepid2::Basis_HGRAD_WEDGE_C2_FEM class.
Wrapper around a Kokkos::View that allows data that is constant or repeating in various logical dimen...
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Line cell.
Header function for Intrepid2::Util class and other utility functions.
Header file for the Intrepid2::Basis_HGRAD_WEDGE_C1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TET_C2_FEM class.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Wedge cell.
Header file for the classes: Intrepid2::RefSubcellParametrization, Intrepid2::RefCellNodes, Intrepid2::RefCellCenter.
Header file for the Intrepid2::Basis_HGRAD_HEX_C1_FEM class.
static bool isSupported(const unsigned cellTopoKey)
Checks if a cell topology has a reference parametrization.
Header file for the Intrepid2::Basis_HGRAD_TET_C1_FEM class.
Contains definitions of custom data types in Intrepid2.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Quadrilateral cell...
Header file for the Intrepid2::Basis_HGRAD_LINE_C1_FEM class.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Tetrahedron cell...
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Tetrahedron cell...
Header file for the Intrepid2::Basis_HGRAD_TET_COMP12_FEM class.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Line cell.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Tetrahedron cell...
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Hexahedron cell...
Header file for the abstract base class Intrepid2::Basis.