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:
 
  157     typedef Kokkos::DynRankView<double,ExecSpaceType> subcellParamViewType;
 
  164     template<
typename outputValueType, 
 
  165              typename pointValueType>
 
  166     static Teuchos::RCP<Basis<ExecSpaceType,outputValueType,pointValueType> >
 
  168       Teuchos::RCP<Basis<ExecSpaceType,outputValueType,pointValueType> > r_val;
 
  170       switch (cellTopo.getKey()) {
 
  189       case shards::Quadrilateral<8>::key: 
 
  190       case shards::Line<3>::key:
 
  191       case shards::Beam<2>::key:
 
  192       case shards::Beam<3>::key:
 
  193       case shards::ShellLine<2>::key:
 
  194       case shards::ShellLine<3>::key:
 
  195       case shards::ShellTriangle<3>::key:
 
  196       case shards::ShellTriangle<6>::key:
 
  197       case shards::ShellQuadrilateral<4>::key:
 
  198       case shards::ShellQuadrilateral<8>::key:
 
  199       case shards::ShellQuadrilateral<9>::key: 
 
  201         INTREPID2_TEST_FOR_EXCEPTION( 
true, std::invalid_argument,
 
  202                                       ">>> ERROR (Intrepid2::CellTools::createHGradBasis): Cell topology not supported.");
 
  212     typedef Kokkos::DynRankView<const double,Kokkos::LayoutRight,Kokkos::HostSpace> referenceNodeDataViewHostType;
 
  214       double line[2][3], line_3[3][3];
 
  215       double triangle[3][3], triangle_4[4][3], triangle_6[6][3];
 
  216       double quadrilateral[4][3], quadrilateral_8[8][3], quadrilateral_9[9][3];
 
  217       double tetrahedron[4][3], tetrahedron_8[8][3], tetrahedron_10[10][3], tetrahedron_11[10][3];
 
  218       double hexahedron[8][3], hexahedron_20[20][3], hexahedron_27[27][3];
 
  219       double pyramid[5][3], pyramid_13[13][3], pyramid_14[14][3];
 
  220       double wedge[6][3], wedge_15[15][3], wedge_18[18][3];
 
  226     typedef Kokkos::DynRankView<double,Kokkos::LayoutRight,ExecSpaceType> referenceNodeDataViewType;
 
  228       referenceNodeDataViewType line, line_3;
 
  229       referenceNodeDataViewType triangle, triangle_4, triangle_6;
 
  230       referenceNodeDataViewType quadrilateral, quadrilateral_8, quadrilateral_9;
 
  231       referenceNodeDataViewType tetrahedron, tetrahedron_8, tetrahedron_10, tetrahedron_11;
 
  232       referenceNodeDataViewType hexahedron, hexahedron_20, hexahedron_27;
 
  233       referenceNodeDataViewType pyramid, pyramid_13, pyramid_14;
 
  234       referenceNodeDataViewType wedge, wedge_15, wedge_18;
 
  240     static bool isReferenceNodeDataSet_;
 
  256       subcellParamViewType dummy;
 
  257       subcellParamViewType lineEdges;  
 
  258       subcellParamViewType triEdges, quadEdges; 
 
  259       subcellParamViewType shellTriEdges, shellQuadEdges; 
 
  260       subcellParamViewType tetEdges, hexEdges, pyrEdges, wedgeEdges; 
 
  261       subcellParamViewType shellTriFaces, shellQuadFaces; 
 
  262       subcellParamViewType tetFaces, hexFaces, pyrFaces, wedgeFaces; 
 
  267     static bool isSubcellParametrizationSet_;
 
  282                                const ordinal_type          subcellDim,
 
  283                                const shards::CellTopology  parentCell );
 
  370     template<
typename jacobianValueType,    
class ...jacobianProperties,
 
  371              typename pointValueType,       
class ...pointProperties,
 
  372              typename worksetCellValueType, 
class ...worksetCellProperties,
 
  373              typename HGradBasisPtrType>
 
  375     setJacobian(       Kokkos::DynRankView<jacobianValueType,jacobianProperties...>       jacobian,
 
  376                  const Kokkos::DynRankView<pointValueType,pointProperties...>             points,
 
  377                  const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
 
  378                  const HGradBasisPtrType basis );
 
  414     template<
typename jacobianValueType,    
class ...jacobianProperties,
 
  415              typename pointValueType,       
class ...pointProperties,
 
  416              typename worksetCellValueType, 
class ...worksetCellProperties>
 
  418     setJacobian(      Kokkos::DynRankView<jacobianValueType,jacobianProperties...> jacobian,
 
  419                  const Kokkos::DynRankView<pointValueType,pointProperties...>       points,
 
  420                  const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
 
  421                  const shards::CellTopology cellTopo ) {
 
  422     auto basis = createHGradBasis<pointValueType,pointValueType>(cellTopo);
 
  439     template<
typename jacobianInvValueType, 
class ...jacobianInvProperties,
 
  440              typename jacobianValueType,    
class ...jacobianProperties>
 
  442     setJacobianInv(       Kokkos::DynRankView<jacobianInvValueType,jacobianInvProperties...> jacobianInv,
 
  443                     const Kokkos::DynRankView<jacobianValueType,jacobianProperties...>       jacobian );
 
  455     template<
typename jacobianDetValueType, 
class ...jacobianDetProperties,
 
  456              typename jacobianValueType,    
class ...jacobianProperties>
 
  458     setJacobianDet(       Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...>  jacobianDet,
 
  459                     const Kokkos::DynRankView<jacobianValueType,jacobianProperties...>        jacobian );
 
  481     template<
typename cellCenterValueType, 
class ...cellCenterProperties,
 
  482              typename cellVertexValueType, 
class ...cellVertexProperties>
 
  485                             Kokkos::DynRankView<cellVertexValueType,cellVertexProperties...> cellVertex,
 
  486                             const shards::CellTopology cell );
 
  498     template<
typename cellVertexValueType, 
class ...cellVertexProperties>
 
  500     getReferenceVertex(       Kokkos::DynRankView<cellVertexValueType,cellVertexProperties...> cellVertex,
 
  501                         const shards::CellTopology cell,
 
  502                         const ordinal_type         vertexOrd );
 
  519     template<
typename subcellVertexValueType, 
class ...subcellVertexProperties>
 
  522                                  const ordinal_type         subcellDim,
 
  523                                  const ordinal_type         subcellOrd,
 
  524                                  const shards::CellTopology parentCell );
 
  543     template<
typename cellNodeValueType, 
class ...cellNodeProperties>
 
  545     getReferenceNode(       Kokkos::DynRankView<cellNodeValueType,cellNodeProperties...> cellNode,
 
  546                       const shards::CellTopology  cell,
 
  547                       const ordinal_type          nodeOrd );
 
  564     template<
typename subcellNodeValueType, 
class ...subcellNodeProperties>
 
  567                               const ordinal_type         subcellDim,
 
  568                               const ordinal_type         subcellOrd,
 
  569                               const shards::CellTopology parentCell );
 
  596     template<
typename refEdgeTangentValueType, 
class ...refEdgeTangentProperties>
 
  598     getReferenceEdgeTangent(       Kokkos::DynRankView<refEdgeTangentValueType,refEdgeTangentProperties...> refEdgeTangent,
 
  599                              const ordinal_type         edgeOrd,
 
  600                              const shards::CellTopology parentCell );
 
  638     template<
typename refFaceTanUValueType, 
class ...refFaceTanUProperties,
 
  639              typename refFaceTanVValueType, 
class ...refFaceTanVProperties>
 
  642                                     Kokkos::DynRankView<refFaceTanVValueType,refFaceTanVProperties...> refFaceTanV,
 
  643                               const ordinal_type         faceOrd,
 
  644                               const shards::CellTopology parentCell );
 
  708     template<
typename refSideNormalValueType, 
class ...refSideNormalProperties>
 
  710     getReferenceSideNormal(       Kokkos::DynRankView<refSideNormalValueType,refSideNormalProperties...> refSideNormal,
 
  711                             const ordinal_type         sideOrd,
 
  712                             const shards::CellTopology parentCell );
 
  752     template<
typename refFaceNormalValueType, 
class ...refFaceNormalProperties>
 
  754     getReferenceFaceNormal(       Kokkos::DynRankView<refFaceNormalValueType,refFaceNormalProperties...> refFaceNormal,
 
  755                             const ordinal_type         faceOrd,
 
  756                             const shards::CellTopology parentCell );
 
  787     template<
typename edgeTangentValueType,     
class ...edgeTangentProperties,
 
  788              typename worksetJacobianValueType, 
class ...worksetJacobianProperties>
 
  791                              const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
 
  792                              const ordinal_type         worksetEdgeOrd,
 
  793                              const shards::CellTopology parentCell );
 
  834     template<
typename faceTanUValueType,        
class ...faceTanUProperties,
 
  835              typename faceTanVValueType,        
class ...faceTanVProperties,
 
  836              typename worksetJacobianValueType, 
class ...worksetJacobianProperties>
 
  839                                    Kokkos::DynRankView<faceTanVValueType,faceTanVProperties...> faceTanV,
 
  840                              const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
 
  841                              const ordinal_type         worksetFaceOrd,
 
  842                              const shards::CellTopology parentCell );
 
  904     template<
typename sideNormalValueType,      
class ...sideNormalProperties,
 
  905              typename worksetJacobianValueType, 
class ...worksetJacobianProperties>
 
  908                             const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
 
  909                             const ordinal_type         worksetSideOrd,
 
  910                             const shards::CellTopology parentCell );
 
  950     template<
typename faceNormalValueType,      
class ...faceNormalProperties,
 
  951              typename worksetJacobianValueType, 
class ...worksetJacobianProperties>
 
  954                             const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
 
  955                             const ordinal_type         worksetFaceOrd,
 
  956                             const shards::CellTopology parentCell );
 
 1000     template<
typename physPointValueType,   
class ...physPointProperties,
 
 1001              typename refPointValueType,    
class ...refPointProperties,
 
 1002              typename worksetCellValueType, 
class ...worksetCellProperties,
 
 1003              typename HGradBasisPtrType>
 
 1005     mapToPhysicalFrame(       Kokkos::DynRankView<physPointValueType,physPointProperties...>     physPoints,
 
 1006                         const Kokkos::DynRankView<refPointValueType,refPointProperties...>       refPoints,
 
 1007                         const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
 
 1008                         const HGradBasisPtrType basis );
 
 1050     template<
typename physPointValueType,   
class ...physPointProperties,
 
 1051              typename refPointValueType,    
class ...refPointProperties,
 
 1052              typename worksetCellValueType, 
class ...worksetCellProperties>
 
 1055                         const Kokkos::DynRankView<refPointValueType,refPointProperties...>       refPoints,
 
 1056                         const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
 
 1057                         const shards::CellTopology cellTopo ) {
 
 1058       auto basis = createHGradBasis<refPointValueType,refPointValueType>(cellTopo);
 
 1082                                const ordinal_type          subcellDim,
 
 1083                                const shards::CellTopology  parentCell );
 
 1136     template<
typename refSubcellPointValueType, 
class ...refSubcellPointProperties,
 
 1137              typename paramPointValueType, 
class ...paramPointProperties>
 
 1139     mapToReferenceSubcell(       Kokkos::DynRankView<refSubcellPointValueType,refSubcellPointProperties...> refSubcellPoints,
 
 1140                            const Kokkos::DynRankView<paramPointValueType,paramPointProperties...>           paramPoints,
 
 1141                            const ordinal_type subcellDim,
 
 1142                            const ordinal_type subcellOrd,
 
 1143                            const shards::CellTopology parentCell );
 
 1156     template<
typename refSubcellPointValueType, 
class ...refSubcellPointProperties,
 
 1157              typename paramPointValueType, 
class ...paramPointProperties>
 
 1159     KOKKOS_INLINE_FUNCTION
 
 1160     mapToReferenceSubcell(       Kokkos::DynRankView<refSubcellPointValueType,refSubcellPointProperties...> refSubcellPoints,
 
 1161                            const Kokkos::DynRankView<paramPointValueType,paramPointProperties...>           paramPoints,
 
 1162                            const subcellParamViewType subcellMap,
 
 1163                            const ordinal_type subcellDim,
 
 1164                            const ordinal_type subcellOrd,
 
 1165                            const ordinal_type parentCellDim);
 
 1218     template<
typename refPointValueType,    
class ...refPointProperties,
 
 1219              typename physPointValueType,   
class ...physPointProperties,
 
 1220              typename worksetCellValueType, 
class ...worksetCellProperties>
 
 1222     mapToReferenceFrame(       Kokkos::DynRankView<refPointValueType,refPointProperties...>    refPoints,
 
 1223                          const Kokkos::DynRankView<physPointValueType,physPointProperties...>  physPoints,
 
 1224                          const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...>   worksetCell,
 
 1225                          const shards::CellTopology cellTopo );
 
 1253     template<
typename refPointValueType,    
class ...refPointProperties,
 
 1254              typename initGuessValueType,   
class ...initGuessProperties,
 
 1255              typename physPointValueType,   
class ...physPointProperties,
 
 1256              typename worksetCellValueType, 
class ...worksetCellProperties,
 
 1257              typename HGradBasisPtrType>
 
 1260                                   const Kokkos::DynRankView<initGuessValueType,initGuessProperties...>     initGuess,
 
 1261                                   const Kokkos::DynRankView<physPointValueType,physPointProperties...>     physPoints,
 
 1262                                   const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
 
 1263                                   const HGradBasisPtrType basis );
 
 1296     template<
typename refPointValueType,    
class ...refPointProperties,
 
 1297              typename initGuessValueType,   
class ...initGuessProperties,
 
 1298              typename physPointValueType,   
class ...physPointProperties,
 
 1299              typename worksetCellValueType, 
class ...worksetCellProperties>
 
 1302                                   const Kokkos::DynRankView<initGuessValueType,initGuessProperties...>     initGuess,
 
 1303                                   const Kokkos::DynRankView<physPointValueType,physPointProperties...>     physPoints,
 
 1304                                   const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
 
 1305                                   const shards::CellTopology cellTopo ) {
 
 1306       auto basis = createHGradBasis<refPointValueType,refPointValueType>(cellTopo);
 
 1394     template<
typename subcvCoordValueType, 
class ...subcvCoordProperties,
 
 1395              typename cellCoordValueType,  
class ...cellCoordProperties>
 
 1397     getSubcvCoords(       Kokkos::DynRankView<subcvCoordValueType,subcvCoordProperties...> subcvCoords, 
 
 1398                     const Kokkos::DynRankView<cellCoordValueType,cellCoordProperties...>   cellCoords,
 
 1399                     const shards::CellTopology primaryCell );
 
 1416     template<
typename pointValueType, 
class ...pointProperties>
 
 1419                          const shards::CellTopology cellTopo,
 
 1420                          const double               thres = threshold() );
 
 1455     template<
typename inCellValueType, 
class ...inCellProperties,                                                       
 
 1456              typename pointValueType, 
class ...pointProperties>
 
 1458                                          const Kokkos::DynRankView<pointValueType,pointProperties...> points,                       
 
 1459                                          const shards::CellTopology cellTopo,                                                       
 
 1460                                          const double thres = threshold() );
 
 1483     template<
typename inCellValueType, 
class ...inCellProperties,                                                       
 
 1484              typename pointValueType, 
class ...pointProperties,                                                        
 
 1485              typename cellWorksetValueType, 
class ...cellWorksetProperties>  
 
 1487                                          const Kokkos::DynRankView<pointValueType,pointProperties...> points,                       
 
 1488                                          const Kokkos::DynRankView<cellWorksetValueType,cellWorksetProperties...> cellWorkset,      
 
 1489                                          const shards::CellTopology cellTopo,                                                       
 
 1490                                          const double thres = threshold() );
 
 1535   template<
typename jacobianViewType,
 
 1536            typename PointViewType,
 
 1537            typename worksetCellViewType>
 
 1539   CellTools_setJacobianArgs( 
const jacobianViewType     jacobian,
 
 1540                              const PointViewType        points,
 
 1541                              const worksetCellViewType  worksetCell,
 
 1542                              const shards::CellTopology cellTopo );
 
 1548   template<
typename jacobianInvViewType,
 
 1549            typename jacobianViewType>
 
 1551   CellTools_setJacobianInvArgs( 
const jacobianInvViewType jacobianInv,
 
 1552                                 const jacobianViewType    jacobian );
 
 1559   template<
typename jacobianDetViewType,
 
 1560            typename jacobianViewType>
 
 1562   CellTools_setJacobianDetArgs( 
const jacobianDetViewType jacobianDet,
 
 1563                                 const jacobianViewType    jacobian );
 
 1572   template<
typename physPointViewType,
 
 1573            typename refPointViewType,
 
 1574            typename worksetCellViewType>
 
 1576   CellTools_mapToPhysicalFrameArgs( 
const physPointViewType    physPoints,
 
 1577                                     const refPointViewType     refPoints,
 
 1578                                     const worksetCellViewType  worksetCell,
 
 1579                                     const shards::CellTopology cellTopo );
 
 1588   template<
typename refPointViewType, 
 
 1589            typename physPointViewType, 
 
 1590            typename worksetCellViewType>
 
 1592   CellTools_mapToReferenceFrameArgs( 
const refPointViewType     refPoints,
 
 1593                                      const physPointViewType    physPoints,
 
 1594                                      const worksetCellViewType  worksetCell,
 
 1595                                      const shards::CellTopology cellTopo );
 
 1606   template<
typename refPointViewType, 
 
 1607            typename initGuessViewType, 
 
 1608            typename physPointViewType, 
 
 1609            typename worksetCellViewType>
 
 1611   CellTools_mapToReferenceFrameInitGuess( 
const refPointViewType     refPoints,
 
 1612                                           const initGuessViewType    initGuess,
 
 1613                                           const physPointViewType    physPoints,
 
 1614                                           const worksetCellViewType  worksetCell,
 
 1615                                           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.