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"
103 template<
typename ExecSpaceType>
115 switch ( cellTopo.getKey() ) {
116 case shards::Line<2>::key:
117 case shards::Line<3>::key:
118 case shards::ShellLine<2>::key:
119 case shards::ShellLine<3>::key:
120 case shards::Beam<2>::key:
121 case shards::Beam<3>::key:
123 case shards::Triangle<3>::key:
125 case shards::Triangle<6>::key:
129 case shards::Quadrilateral<4>::key:
130 case shards::Quadrilateral<8>::key:
131 case shards::Quadrilateral<9>::key:
136 case shards::Tetrahedron<4>::key:
138 case shards::Tetrahedron<10>::key:
141 case shards::Hexahedron<8>::key:
142 case shards::Hexahedron<20>::key:
143 case shards::Hexahedron<27>::key:
145 case shards::Pyramid<5>::key:
149 case shards::Wedge<6>::key:
151 case shards::Wedge<18>::key:
162 template<
typename outputValueType,
163 typename pointValueType>
164 static Teuchos::RCP<Basis<ExecSpaceType,outputValueType,pointValueType> >
166 Teuchos::RCP<Basis<ExecSpaceType,outputValueType,pointValueType> > r_val;
168 switch (cellTopo.getKey()) {
187 case shards::Quadrilateral<8>::key:
188 case shards::Line<3>::key:
189 case shards::Beam<2>::key:
190 case shards::Beam<3>::key:
191 case shards::ShellLine<2>::key:
192 case shards::ShellLine<3>::key:
193 case shards::ShellTriangle<3>::key:
194 case shards::ShellTriangle<6>::key:
195 case shards::ShellQuadrilateral<4>::key:
196 case shards::ShellQuadrilateral<8>::key:
197 case shards::ShellQuadrilateral<9>::key:
199 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
200 ">>> ERROR (Intrepid2::CellTools::createHGradBasis): Cell topology not supported.");
210 typedef Kokkos::DynRankView<const double,Kokkos::LayoutRight,Kokkos::HostSpace> referenceNodeDataViewHostType;
212 double line[2][3], line_3[3][3];
213 double triangle[3][3], triangle_4[4][3], triangle_6[6][3];
214 double quadrilateral[4][3], quadrilateral_8[8][3], quadrilateral_9[9][3];
215 double tetrahedron[4][3], tetrahedron_8[8][3], tetrahedron_10[10][3], tetrahedron_11[10][3];
216 double hexahedron[8][3], hexahedron_20[20][3], hexahedron_27[27][3];
217 double pyramid[5][3], pyramid_13[13][3], pyramid_14[14][3];
218 double wedge[6][3], wedge_15[15][3], wedge_18[18][3];
224 typedef Kokkos::DynRankView<double,Kokkos::LayoutRight,ExecSpaceType> referenceNodeDataViewType;
226 referenceNodeDataViewType line, line_3;
227 referenceNodeDataViewType triangle, triangle_4, triangle_6;
228 referenceNodeDataViewType quadrilateral, quadrilateral_8, quadrilateral_9;
229 referenceNodeDataViewType tetrahedron, tetrahedron_8, tetrahedron_10, tetrahedron_11;
230 referenceNodeDataViewType hexahedron, hexahedron_20, hexahedron_27;
231 referenceNodeDataViewType pyramid, pyramid_13, pyramid_14;
232 referenceNodeDataViewType wedge, wedge_15, wedge_18;
238 static bool isReferenceNodeDataSet_;
253 typedef Kokkos::DynRankView<double,ExecSpaceType> subcellParamViewType;
255 subcellParamViewType dummy;
256 subcellParamViewType lineEdges;
257 subcellParamViewType triEdges, quadEdges;
258 subcellParamViewType shellTriEdges, shellQuadEdges;
259 subcellParamViewType tetEdges, hexEdges, pyrEdges, wedgeEdges;
260 subcellParamViewType shellTriFaces, shellQuadFaces;
261 subcellParamViewType tetFaces, hexFaces, pyrFaces, wedgeFaces;
266 static bool isSubcellParametrizationSet_;
316 const ordinal_type subcellDim,
317 const shards::CellTopology parentCell );
368 template<
typename jacobianValueType,
class ...jacobianProperties,
369 typename pointValueType,
class ...pointProperties,
370 typename worksetCellValueType,
class ...worksetCellProperties,
371 typename HGradBasisPtrType>
373 setJacobian( Kokkos::DynRankView<jacobianValueType,jacobianProperties...> jacobian,
374 const Kokkos::DynRankView<pointValueType,pointProperties...> points,
375 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
376 const HGradBasisPtrType basis );
412 template<
typename jacobianValueType,
class ...jacobianProperties,
413 typename pointValueType,
class ...pointProperties,
414 typename worksetCellValueType,
class ...worksetCellProperties>
416 setJacobian( Kokkos::DynRankView<jacobianValueType,jacobianProperties...> jacobian,
417 const Kokkos::DynRankView<pointValueType,pointProperties...> points,
418 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
419 const shards::CellTopology cellTopo ) {
420 auto basis = createHGradBasis<pointValueType,pointValueType>(cellTopo);
437 template<
typename jacobianInvValueType,
class ...jacobianInvProperties,
438 typename jacobianValueType,
class ...jacobianProperties>
440 setJacobianInv( Kokkos::DynRankView<jacobianInvValueType,jacobianInvProperties...> jacobianInv,
441 const Kokkos::DynRankView<jacobianValueType,jacobianProperties...> jacobian );
453 template<
typename jacobianDetValueType,
class ...jacobianDetProperties,
454 typename jacobianValueType,
class ...jacobianProperties>
456 setJacobianDet( Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
457 const Kokkos::DynRankView<jacobianValueType,jacobianProperties...> jacobian );
479 template<
typename cellCenterValueType,
class ...cellCenterProperties,
480 typename cellVertexValueType,
class ...cellVertexProperties>
483 Kokkos::DynRankView<cellVertexValueType,cellVertexProperties...> cellVertex,
484 const shards::CellTopology cell );
496 template<
typename cellVertexValueType,
class ...cellVertexProperties>
498 getReferenceVertex( Kokkos::DynRankView<cellVertexValueType,cellVertexProperties...> cellVertex,
499 const shards::CellTopology cell,
500 const ordinal_type vertexOrd );
517 template<
typename subcellVertexValueType,
class ...subcellVertexProperties>
520 const ordinal_type subcellDim,
521 const ordinal_type subcellOrd,
522 const shards::CellTopology parentCell );
541 template<
typename cellNodeValueType,
class ...cellNodeProperties>
543 getReferenceNode( Kokkos::DynRankView<cellNodeValueType,cellNodeProperties...> cellNode,
544 const shards::CellTopology cell,
545 const ordinal_type nodeOrd );
562 template<
typename subcellNodeValueType,
class ...subcellNodeProperties>
565 const ordinal_type subcellDim,
566 const ordinal_type subcellOrd,
567 const shards::CellTopology parentCell );
594 template<
typename refEdgeTangentValueType,
class ...refEdgeTangentProperties>
596 getReferenceEdgeTangent( Kokkos::DynRankView<refEdgeTangentValueType,refEdgeTangentProperties...> refEdgeTangent,
597 const ordinal_type edgeOrd,
598 const shards::CellTopology parentCell );
636 template<
typename refFaceTanUValueType,
class ...refFaceTanUProperties,
637 typename refFaceTanVValueType,
class ...refFaceTanVProperties>
640 Kokkos::DynRankView<refFaceTanVValueType,refFaceTanVProperties...> refFaceTanV,
641 const ordinal_type faceOrd,
642 const shards::CellTopology parentCell );
706 template<
typename refSideNormalValueType,
class ...refSideNormalProperties>
708 getReferenceSideNormal( Kokkos::DynRankView<refSideNormalValueType,refSideNormalProperties...> refSideNormal,
709 const ordinal_type sideOrd,
710 const shards::CellTopology parentCell );
750 template<
typename refFaceNormalValueType,
class ...refFaceNormalProperties>
752 getReferenceFaceNormal( Kokkos::DynRankView<refFaceNormalValueType,refFaceNormalProperties...> refFaceNormal,
753 const ordinal_type faceOrd,
754 const shards::CellTopology parentCell );
785 template<
typename edgeTangentValueType,
class ...edgeTangentProperties,
786 typename worksetJacobianValueType,
class ...worksetJacobianProperties>
789 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
790 const ordinal_type worksetEdgeOrd,
791 const shards::CellTopology parentCell );
832 template<
typename faceTanUValueType,
class ...faceTanUProperties,
833 typename faceTanVValueType,
class ...faceTanVProperties,
834 typename worksetJacobianValueType,
class ...worksetJacobianProperties>
837 Kokkos::DynRankView<faceTanVValueType,faceTanVProperties...> faceTanV,
838 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
839 const ordinal_type worksetFaceOrd,
840 const shards::CellTopology parentCell );
902 template<
typename sideNormalValueType,
class ...sideNormalProperties,
903 typename worksetJacobianValueType,
class ...worksetJacobianProperties>
906 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
907 const ordinal_type worksetSideOrd,
908 const shards::CellTopology parentCell );
948 template<
typename faceNormalValueType,
class ...faceNormalProperties,
949 typename worksetJacobianValueType,
class ...worksetJacobianProperties>
952 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
953 const ordinal_type worksetFaceOrd,
954 const shards::CellTopology parentCell );
998 template<
typename physPointValueType,
class ...physPointProperties,
999 typename refPointValueType,
class ...refPointProperties,
1000 typename worksetCellValueType,
class ...worksetCellProperties,
1001 typename HGradBasisPtrType>
1003 mapToPhysicalFrame( Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1004 const Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
1005 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1006 const HGradBasisPtrType basis );
1048 template<
typename physPointValueType,
class ...physPointProperties,
1049 typename refPointValueType,
class ...refPointProperties,
1050 typename worksetCellValueType,
class ...worksetCellProperties>
1053 const Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
1054 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1055 const shards::CellTopology cellTopo ) {
1056 auto basis = createHGradBasis<refPointValueType,refPointValueType>(cellTopo);
1080 const ordinal_type subcellDim,
1081 const shards::CellTopology parentCell );
1134 template<
typename refSubcellPointValueType,
class ...refSubcellPointProperties,
1135 typename paramPointValueType,
class ...paramPointProperties>
1137 mapToReferenceSubcell( Kokkos::DynRankView<refSubcellPointValueType,refSubcellPointProperties...> refSubcellPoints,
1138 const Kokkos::DynRankView<paramPointValueType,paramPointProperties...> paramPoints,
1139 const ordinal_type subcellDim,
1140 const ordinal_type subcellOrd,
1141 const shards::CellTopology parentCell );
1194 template<
typename refPointValueType,
class ...refPointProperties,
1195 typename physPointValueType,
class ...physPointProperties,
1196 typename worksetCellValueType,
class ...worksetCellProperties>
1198 mapToReferenceFrame( Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
1199 const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1200 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1201 const shards::CellTopology cellTopo );
1229 template<
typename refPointValueType,
class ...refPointProperties,
1230 typename initGuessValueType,
class ...initGuessProperties,
1231 typename physPointValueType,
class ...physPointProperties,
1232 typename worksetCellValueType,
class ...worksetCellProperties,
1233 typename HGradBasisPtrType>
1236 const Kokkos::DynRankView<initGuessValueType,initGuessProperties...> initGuess,
1237 const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1238 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1239 const HGradBasisPtrType basis );
1272 template<
typename refPointValueType,
class ...refPointProperties,
1273 typename initGuessValueType,
class ...initGuessProperties,
1274 typename physPointValueType,
class ...physPointProperties,
1275 typename worksetCellValueType,
class ...worksetCellProperties>
1278 const Kokkos::DynRankView<initGuessValueType,initGuessProperties...> initGuess,
1279 const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1280 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1281 const shards::CellTopology cellTopo ) {
1282 auto basis = createHGradBasis<refPointValueType,refPointValueType>(cellTopo);
1370 template<
typename subcvCoordValueType,
class ...subcvCoordProperties,
1371 typename cellCoordValueType,
class ...cellCoordProperties>
1373 getSubcvCoords( Kokkos::DynRankView<subcvCoordValueType,subcvCoordProperties...> subcvCoords,
1374 const Kokkos::DynRankView<cellCoordValueType,cellCoordProperties...> cellCoords,
1375 const shards::CellTopology primaryCell );
1392 template<
typename pointValueType,
class ...pointProperties>
1395 const shards::CellTopology cellTopo,
1396 const double thres = threshold() );
1431 template<
typename inCellValueType,
class ...inCellProperties,
1432 typename pointValueType,
class ...pointProperties>
1434 const Kokkos::DynRankView<pointValueType,pointProperties...> points,
1435 const shards::CellTopology cellTopo,
1436 const double thres = threshold() );
1459 template<
typename inCellValueType,
class ...inCellProperties,
1460 typename pointValueType,
class ...pointProperties,
1461 typename cellWorksetValueType,
class ...cellWorksetProperties>
1463 const Kokkos::DynRankView<pointValueType,pointProperties...> points,
1464 const Kokkos::DynRankView<cellWorksetValueType,cellWorksetProperties...> cellWorkset,
1465 const shards::CellTopology cellTopo,
1466 const double thres = threshold() );
1511 template<
typename jacobianViewType,
1512 typename PointViewType,
1513 typename worksetCellViewType>
1515 CellTools_setJacobianArgs(
const jacobianViewType jacobian,
1516 const PointViewType points,
1517 const worksetCellViewType worksetCell,
1518 const shards::CellTopology cellTopo );
1524 template<
typename jacobianInvViewType,
1525 typename jacobianViewType>
1527 CellTools_setJacobianInvArgs(
const jacobianInvViewType jacobianInv,
1528 const jacobianViewType jacobian );
1535 template<
typename jacobianDetViewType,
1536 typename jacobianViewType>
1538 CellTools_setJacobianDetArgs(
const jacobianDetViewType jacobianDet,
1539 const jacobianViewType jacobian );
1548 template<
typename physPointViewType,
1549 typename refPointViewType,
1550 typename worksetCellViewType>
1552 CellTools_mapToPhysicalFrameArgs(
const physPointViewType physPoints,
1553 const refPointViewType refPoints,
1554 const worksetCellViewType worksetCell,
1555 const shards::CellTopology cellTopo );
1564 template<
typename refPointViewType,
1565 typename physPointViewType,
1566 typename worksetCellViewType>
1568 CellTools_mapToReferenceFrameArgs(
const refPointViewType refPoints,
1569 const physPointViewType physPoints,
1570 const worksetCellViewType worksetCell,
1571 const shards::CellTopology cellTopo );
1582 template<
typename refPointViewType,
1583 typename initGuessViewType,
1584 typename physPointViewType,
1585 typename worksetCellViewType>
1587 CellTools_mapToReferenceFrameInitGuess(
const refPointViewType refPoints,
1588 const initGuessViewType initGuess,
1589 const physPointViewType physPoints,
1590 const worksetCellViewType worksetCell,
1591 const shards::CellTopology cellTopo );
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Wedge cell.
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.
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 2 on Hexahedron cell...
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...
Header file for the Intrepid2::Basis_HGRAD_HEX_C2_FEM class.
Header file for the Intrepid2::Basis_HGRAD_PYR_C1_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.
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.
Header file for the Intrepid2::Basis_HGRAD_HEX_C1_FEM class.
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.