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_;
318 const ordinal_type subcellDim,
319 const shards::CellTopology parentCell );
334 const ordinal_type subcellDim,
335 const shards::CellTopology parentCell );
386 template<
typename jacobianValueType,
class ...jacobianProperties,
387 typename pointValueType,
class ...pointProperties,
388 typename worksetCellValueType,
class ...worksetCellProperties,
389 typename HGradBasisPtrType>
391 setJacobian( Kokkos::DynRankView<jacobianValueType,jacobianProperties...> jacobian,
392 const Kokkos::DynRankView<pointValueType,pointProperties...> points,
393 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
394 const HGradBasisPtrType basis );
430 template<
typename jacobianValueType,
class ...jacobianProperties,
431 typename pointValueType,
class ...pointProperties,
432 typename worksetCellValueType,
class ...worksetCellProperties>
434 setJacobian( Kokkos::DynRankView<jacobianValueType,jacobianProperties...> jacobian,
435 const Kokkos::DynRankView<pointValueType,pointProperties...> points,
436 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
437 const shards::CellTopology cellTopo ) {
438 auto basis = createHGradBasis<pointValueType,pointValueType>(cellTopo);
455 template<
typename jacobianInvValueType,
class ...jacobianInvProperties,
456 typename jacobianValueType,
class ...jacobianProperties>
458 setJacobianInv( Kokkos::DynRankView<jacobianInvValueType,jacobianInvProperties...> jacobianInv,
459 const Kokkos::DynRankView<jacobianValueType,jacobianProperties...> jacobian );
471 template<
typename jacobianDetValueType,
class ...jacobianDetProperties,
472 typename jacobianValueType,
class ...jacobianProperties>
474 setJacobianDet( Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
475 const Kokkos::DynRankView<jacobianValueType,jacobianProperties...> jacobian );
497 template<
typename cellCenterValueType,
class ...cellCenterProperties,
498 typename cellVertexValueType,
class ...cellVertexProperties>
501 Kokkos::DynRankView<cellVertexValueType,cellVertexProperties...> cellVertex,
502 const shards::CellTopology cell );
514 template<
typename cellVertexValueType,
class ...cellVertexProperties>
516 getReferenceVertex( Kokkos::DynRankView<cellVertexValueType,cellVertexProperties...> cellVertex,
517 const shards::CellTopology cell,
518 const ordinal_type vertexOrd );
535 template<
typename subcellVertexValueType,
class ...subcellVertexProperties>
538 const ordinal_type subcellDim,
539 const ordinal_type subcellOrd,
540 const shards::CellTopology parentCell );
559 template<
typename cellNodeValueType,
class ...cellNodeProperties>
561 getReferenceNode( Kokkos::DynRankView<cellNodeValueType,cellNodeProperties...> cellNode,
562 const shards::CellTopology cell,
563 const ordinal_type nodeOrd );
580 template<
typename subcellNodeValueType,
class ...subcellNodeProperties>
583 const ordinal_type subcellDim,
584 const ordinal_type subcellOrd,
585 const shards::CellTopology parentCell );
612 template<
typename refEdgeTangentValueType,
class ...refEdgeTangentProperties>
614 getReferenceEdgeTangent( Kokkos::DynRankView<refEdgeTangentValueType,refEdgeTangentProperties...> refEdgeTangent,
615 const ordinal_type edgeOrd,
616 const shards::CellTopology parentCell );
654 template<
typename refFaceTanUValueType,
class ...refFaceTanUProperties,
655 typename refFaceTanVValueType,
class ...refFaceTanVProperties>
658 Kokkos::DynRankView<refFaceTanVValueType,refFaceTanVProperties...> refFaceTanV,
659 const ordinal_type faceOrd,
660 const shards::CellTopology parentCell );
724 template<
typename refSideNormalValueType,
class ...refSideNormalProperties>
726 getReferenceSideNormal( Kokkos::DynRankView<refSideNormalValueType,refSideNormalProperties...> refSideNormal,
727 const ordinal_type sideOrd,
728 const shards::CellTopology parentCell );
768 template<
typename refFaceNormalValueType,
class ...refFaceNormalProperties>
770 getReferenceFaceNormal( Kokkos::DynRankView<refFaceNormalValueType,refFaceNormalProperties...> refFaceNormal,
771 const ordinal_type faceOrd,
772 const shards::CellTopology parentCell );
803 template<
typename edgeTangentValueType,
class ...edgeTangentProperties,
804 typename worksetJacobianValueType,
class ...worksetJacobianProperties>
807 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
808 const ordinal_type worksetEdgeOrd,
809 const shards::CellTopology parentCell );
850 template<
typename faceTanUValueType,
class ...faceTanUProperties,
851 typename faceTanVValueType,
class ...faceTanVProperties,
852 typename worksetJacobianValueType,
class ...worksetJacobianProperties>
855 Kokkos::DynRankView<faceTanVValueType,faceTanVProperties...> faceTanV,
856 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
857 const ordinal_type worksetFaceOrd,
858 const shards::CellTopology parentCell );
920 template<
typename sideNormalValueType,
class ...sideNormalProperties,
921 typename worksetJacobianValueType,
class ...worksetJacobianProperties>
924 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
925 const ordinal_type worksetSideOrd,
926 const shards::CellTopology parentCell );
966 template<
typename faceNormalValueType,
class ...faceNormalProperties,
967 typename worksetJacobianValueType,
class ...worksetJacobianProperties>
970 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
971 const ordinal_type worksetFaceOrd,
972 const shards::CellTopology parentCell );
1016 template<
typename physPointValueType,
class ...physPointProperties,
1017 typename refPointValueType,
class ...refPointProperties,
1018 typename worksetCellValueType,
class ...worksetCellProperties,
1019 typename HGradBasisPtrType>
1021 mapToPhysicalFrame( Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1022 const Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
1023 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1024 const HGradBasisPtrType basis );
1066 template<
typename physPointValueType,
class ...physPointProperties,
1067 typename refPointValueType,
class ...refPointProperties,
1068 typename worksetCellValueType,
class ...worksetCellProperties>
1071 const Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
1072 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1073 const shards::CellTopology cellTopo ) {
1074 auto basis = createHGradBasis<refPointValueType,refPointValueType>(cellTopo);
1131 template<
typename refSubcellPointValueType,
class ...refSubcellPointProperties,
1132 typename paramPointValueType,
class ...paramPointProperties>
1134 mapToReferenceSubcell( Kokkos::DynRankView<refSubcellPointValueType,refSubcellPointProperties...> refSubcellPoints,
1135 const Kokkos::DynRankView<paramPointValueType,paramPointProperties...> paramPoints,
1136 const ordinal_type subcellDim,
1137 const ordinal_type subcellOrd,
1138 const shards::CellTopology parentCell );
1191 template<
typename refPointValueType,
class ...refPointProperties,
1192 typename physPointValueType,
class ...physPointProperties,
1193 typename worksetCellValueType,
class ...worksetCellProperties>
1195 mapToReferenceFrame( Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
1196 const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1197 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1198 const shards::CellTopology cellTopo );
1226 template<
typename refPointValueType,
class ...refPointProperties,
1227 typename initGuessValueType,
class ...initGuessProperties,
1228 typename physPointValueType,
class ...physPointProperties,
1229 typename worksetCellValueType,
class ...worksetCellProperties,
1230 typename HGradBasisPtrType>
1233 const Kokkos::DynRankView<initGuessValueType,initGuessProperties...> initGuess,
1234 const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1235 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1236 const HGradBasisPtrType basis );
1269 template<
typename refPointValueType,
class ...refPointProperties,
1270 typename initGuessValueType,
class ...initGuessProperties,
1271 typename physPointValueType,
class ...physPointProperties,
1272 typename worksetCellValueType,
class ...worksetCellProperties>
1275 const Kokkos::DynRankView<initGuessValueType,initGuessProperties...> initGuess,
1276 const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1277 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1278 const shards::CellTopology cellTopo ) {
1279 auto basis = createHGradBasis<refPointValueType,refPointValueType>(cellTopo);
1367 template<
typename subcvCoordValueType,
class ...subcvCoordProperties,
1368 typename cellCoordValueType,
class ...cellCoordProperties>
1370 getSubcvCoords( Kokkos::DynRankView<subcvCoordValueType,subcvCoordProperties...> subcvCoords,
1371 const Kokkos::DynRankView<cellCoordValueType,cellCoordProperties...> cellCoords,
1372 const shards::CellTopology primaryCell );
1389 template<
typename pointValueType,
class ...pointProperties>
1392 const shards::CellTopology cellTopo,
1393 const double thres = threshold() );
1428 template<
typename inCellValueType,
class ...inCellProperties,
1429 typename pointValueType,
class ...pointProperties>
1431 const Kokkos::DynRankView<pointValueType,pointProperties...> points,
1432 const shards::CellTopology cellTopo,
1433 const double thres = threshold() );
1456 template<
typename inCellValueType,
class ...inCellProperties,
1457 typename pointValueType,
class ...pointProperties,
1458 typename cellWorksetValueType,
class ...cellWorksetProperties>
1460 const Kokkos::DynRankView<pointValueType,pointProperties...> points,
1461 const Kokkos::DynRankView<cellWorksetValueType,cellWorksetProperties...> cellWorkset,
1462 const shards::CellTopology cellTopo,
1463 const double thres = threshold() );
1508 template<
typename jacobianViewType,
1509 typename pointViewType,
1510 typename worksetCellViewType>
1512 CellTools_setJacobianArgs(
const jacobianViewType jacobian,
1513 const pointViewType points,
1514 const worksetCellViewType worksetCell,
1515 const shards::CellTopology cellTopo );
1521 template<
typename jacobianInvViewType,
1522 typename jacobianViewType>
1524 CellTools_setJacobianInvArgs(
const jacobianInvViewType jacobianInv,
1525 const jacobianViewType jacobian );
1532 template<
typename jacobianDetViewType,
1533 typename jacobianViewType>
1535 CellTools_setJacobianDetArgs(
const jacobianDetViewType jacobianDet,
1536 const jacobianViewType jacobian );
1545 template<
typename physPointViewType,
1546 typename refPointViewType,
1547 typename worksetCellViewType>
1549 CellTools_mapToPhysicalFrameArgs(
const physPointViewType physPoints,
1550 const refPointViewType refPoints,
1551 const worksetCellViewType worksetCell,
1552 const shards::CellTopology cellTopo );
1561 template<
typename refPointViewType,
1562 typename physPointViewType,
1563 typename worksetCellViewType>
1565 CellTools_mapToReferenceFrameArgs(
const refPointViewType refPoints,
1566 const physPointViewType physPoints,
1567 const worksetCellViewType worksetCell,
1568 const shards::CellTopology cellTopo );
1579 template<
typename refPointViewType,
1580 typename initGuessViewType,
1581 typename physPointViewType,
1582 typename worksetCellViewType>
1584 CellTools_mapToReferenceFrameInitGuess(
const refPointViewType refPoints,
1585 const initGuessViewType initGuess,
1586 const physPointViewType physPoints,
1587 const worksetCellViewType worksetCell,
1588 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.