Intrepid2
Intrepid2_CellTools.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Intrepid2 Package
4 //
5 // Copyright 2007 NTESS and the Intrepid2 contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
16 #ifndef __INTREPID2_CELLTOOLS_HPP__
17 #define __INTREPID2_CELLTOOLS_HPP__
18 
19 #include "Intrepid2_ConfigDefs.hpp"
20 
21 #include "Shards_CellTopology.hpp"
22 #include "Shards_BasicTopologies.hpp"
23 
24 #include "Teuchos_RCP.hpp"
25 
26 #include "Intrepid2_Types.hpp"
27 #include "Intrepid2_Utils.hpp"
28 
30 
31 #include "Intrepid2_Basis.hpp"
32 
35 
38 
41 
45 
48 
51 
54 
55 #include "Intrepid2_Data.hpp"
56 #include "Intrepid2_CellData.hpp"
57 
58 namespace Intrepid2 {
59 
60  //============================================================================================//
61  // //
62  // CellTools //
63  // //
64  //============================================================================================//
65 
76  template<typename DeviceType>
77  class CellTools {
78  using ExecSpaceType = typename DeviceType::execution_space;
79  using MemSpaceType = typename DeviceType::memory_space;
80  public:
81 
86  inline
87  static bool
88  hasReferenceCell( const shards::CellTopology cellTopo ) {
90  }
91 
92  private:
93 
97  template<typename outputValueType,
98  typename pointValueType>
99  static Teuchos::RCP<Basis<DeviceType,outputValueType,pointValueType> >
100  createHGradBasis( const shards::CellTopology cellTopo ) {
101  Teuchos::RCP<Basis<DeviceType,outputValueType,pointValueType> > r_val;
102 
103  switch (cellTopo.getKey()) {
104  case shards::Line<2>::key: r_val = Teuchos::rcp(new Basis_HGRAD_LINE_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
105  case shards::Line<3>::key: r_val = Teuchos::rcp(new Basis_HGRAD_LINE_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
106  case shards::Triangle<3>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TRI_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
107  case shards::Quadrilateral<4>::key: r_val = Teuchos::rcp(new Basis_HGRAD_QUAD_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
108  case shards::Tetrahedron<4>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TET_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
109  case shards::Hexahedron<8>::key: r_val = Teuchos::rcp(new Basis_HGRAD_HEX_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
110  case shards::Wedge<6>::key: r_val = Teuchos::rcp(new Basis_HGRAD_WEDGE_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
111  case shards::Pyramid<5>::key: r_val = Teuchos::rcp(new Basis_HGRAD_PYR_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
112 
113  case shards::Triangle<6>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TRI_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
114  case shards::Quadrilateral<8>::key: r_val = Teuchos::rcp(new Basis_HGRAD_QUAD_I2_FEM <DeviceType,outputValueType,pointValueType>()); break;
115  case shards::Quadrilateral<9>::key: r_val = Teuchos::rcp(new Basis_HGRAD_QUAD_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
116  case shards::Tetrahedron<10>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TET_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
117  case shards::Tetrahedron<11>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TET_COMP12_FEM<DeviceType,outputValueType,pointValueType>()); break;
118  case shards::Hexahedron<20>::key: r_val = Teuchos::rcp(new Basis_HGRAD_HEX_I2_FEM <DeviceType,outputValueType,pointValueType>()); break;
119  case shards::Hexahedron<27>::key: r_val = Teuchos::rcp(new Basis_HGRAD_HEX_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
120  case shards::Wedge<15>::key: r_val = Teuchos::rcp(new Basis_HGRAD_WEDGE_I2_FEM <DeviceType,outputValueType,pointValueType>()); break;
121  case shards::Wedge<18>::key: r_val = Teuchos::rcp(new Basis_HGRAD_WEDGE_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
122  case shards::Pyramid<13>::key: r_val = Teuchos::rcp(new Basis_HGRAD_PYR_I2_FEM <DeviceType,outputValueType,pointValueType>()); break;
123 
124  case shards::Beam<2>::key:
125  case shards::Beam<3>::key:
126  case shards::ShellLine<2>::key:
127  case shards::ShellLine<3>::key:
128  case shards::ShellTriangle<3>::key:
129  case shards::ShellTriangle<6>::key:
130  case shards::ShellQuadrilateral<4>::key:
131  case shards::ShellQuadrilateral<8>::key:
132  case shards::ShellQuadrilateral<9>::key:
133  default: {
134  INTREPID2_TEST_FOR_EXCEPTION( true, std::invalid_argument,
135  ">>> ERROR (Intrepid2::CellTools::createHGradBasis): Cell topology not supported.");
136  }
137  }
138  return r_val;
139  }
140 
141 public:
142 
145  CellTools() = default;
146 
149  ~CellTools() = default;
150 
151  //============================================================================================//
152  // //
153  // Jacobian, inverse Jacobian and Jacobian determinant //
154  // //
155  //============================================================================================//
156 
192  template<typename JacobianViewType,
193  typename PointViewType,
194  typename WorksetType,
195  typename HGradBasisType>
196  static void
197  setJacobian( JacobianViewType jacobian,
198  const PointViewType points,
199  const WorksetType worksetCell,
200  const Teuchos::RCP<HGradBasisType> basis,
201  const int startCell=0, const int endCell=-1);
202 
237  template<typename JacobianViewType,
238  typename BasisGradientsType,
239  typename WorksetType>
240  static void
241  setJacobian( JacobianViewType jacobian,
242  const WorksetType worksetCell,
243  const BasisGradientsType gradients,
244  const int startCell=0, const int endCell=-1);
245 
280  template<typename JacobianViewType,
281  typename PointViewType,
282  typename WorksetCellViewType>
283  static void
284  setJacobian( JacobianViewType jacobian,
285  const PointViewType points,
286  const WorksetCellViewType worksetCell,
287  const shards::CellTopology cellTopo ) {
288  using nonConstPointValueType = typename PointViewType::non_const_value_type;
289  auto basis = createHGradBasis<nonConstPointValueType,nonConstPointValueType>(cellTopo);
290  setJacobian(jacobian,
291  points,
292  worksetCell,
293  basis);
294  }
295 
306  template<typename JacobianInvViewType,
307  typename JacobianViewType>
308  static void
309  setJacobianInv( JacobianInvViewType jacobianInv,
310  const JacobianViewType jacobian );
311 
322  template<typename JacobianDetViewType,
323  typename JacobianViewType>
324  static void
325  setJacobianDet( JacobianDetViewType jacobianDet,
326  const JacobianViewType jacobian );
327 
333  template<class PointScalar>
335 
341  template<class PointScalar>
343 
349  template<class PointScalar>
350  static void setJacobianDet( Data<PointScalar,DeviceType> & jacobianDet,
351  const Data<PointScalar,DeviceType> & jacobian);
352 
358  template<class PointScalar>
359  static void setJacobianDetInv( Data<PointScalar,DeviceType> & jacobianDetInv,
360  const Data<PointScalar,DeviceType> & jacobian);
361 
367  template<class PointScalar>
368  static void setJacobianInv( Data<PointScalar,DeviceType> & jacobianInv,
369  const Data<PointScalar,DeviceType> & jacobian);
370 
377  template<class PointScalar>
378  static void setJacobianDividedByDet( Data<PointScalar,DeviceType> & jacobianDividedByDet,
379  const Data<PointScalar,DeviceType> & jacobian,
380  const Data<PointScalar,DeviceType> & jacobianDetInv);
381 
382  //============================================================================================//
383  // //
384  // Node information //
385  // //
386  //============================================================================================//
387 
388  // the node information can be used inside of kokkos functor and needs kokkos inline and
389  // exception should be an abort. for now, let's not decorate
390 
398  template<typename cellCenterValueType, class ...cellCenterProperties>
399  static void
400  getReferenceCellCenter( Kokkos::DynRankView<cellCenterValueType,cellCenterProperties...> cellCenter,
401  const shards::CellTopology cell );
402 
411  template<typename cellVertexValueType, class ...cellVertexProperties>
412  static void
413  getReferenceVertex( Kokkos::DynRankView<cellVertexValueType,cellVertexProperties...> cellVertex,
414  const shards::CellTopology cell,
415  const ordinal_type vertexOrd );
416 
417 
432  template<typename subcellVertexValueType, class ...subcellVertexProperties>
433  static void
434  getReferenceSubcellVertices( Kokkos::DynRankView<subcellVertexValueType,subcellVertexProperties...> subcellVertices,
435  const ordinal_type subcellDim,
436  const ordinal_type subcellOrd,
437  const shards::CellTopology parentCell );
438 
439 
440 
456  template<typename cellNodeValueType, class ...cellNodeProperties>
457  static void
458  getReferenceNode( Kokkos::DynRankView<cellNodeValueType,cellNodeProperties...> cellNode,
459  const shards::CellTopology cell,
460  const ordinal_type nodeOrd );
461 
462 
463 
477  template<typename SubcellNodeViewType>
478  static void
479  getReferenceSubcellNodes( SubcellNodeViewType subcellNodes,
480  const ordinal_type subcellDim,
481  const ordinal_type subcellOrd,
482  const shards::CellTopology parentCell );
483 
509  template<typename RefEdgeTangentViewType>
510  static void
511  getReferenceEdgeTangent( RefEdgeTangentViewType refEdgeTangent,
512  const ordinal_type edgeOrd,
513  const shards::CellTopology parentCell );
514 
551  template<typename RefFaceTanViewType>
552  static void
553  getReferenceFaceTangents( RefFaceTanViewType refFaceTanU,
554  RefFaceTanViewType refFaceTanV,
555  const ordinal_type faceOrd,
556  const shards::CellTopology parentCell );
557 
620  template<typename RefSideNormalViewType>
621  static void
622  getReferenceSideNormal( RefSideNormalViewType refSideNormal,
623  const ordinal_type sideOrd,
624  const shards::CellTopology parentCell );
625 
664  template<typename RefFaceNormalViewType>
665  static void
666  getReferenceFaceNormal( RefFaceNormalViewType refFaceNormal,
667  const ordinal_type faceOrd,
668  const shards::CellTopology parentCell );
669 
699  template<typename edgeTangentValueType, class ...edgeTangentProperties,
700  typename worksetJacobianValueType, class ...worksetJacobianProperties>
701  static void
702  getPhysicalEdgeTangents( Kokkos::DynRankView<edgeTangentValueType,edgeTangentProperties...> edgeTangents,
703  const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
704  const ordinal_type worksetEdgeOrd,
705  const shards::CellTopology parentCell );
706 
707 
719  template<typename edgeTangentValueType, class ...edgeTangentProperties,
720  typename worksetJacobianValueType, class ...worksetJacobianProperties,
721  typename edgeOrdValueType, class ...edgeOrdProperties>
722  static void
723  getPhysicalEdgeTangents( Kokkos::DynRankView<edgeTangentValueType,edgeTangentProperties...> edgeTangents,
724  const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
725  const Kokkos::DynRankView<edgeOrdValueType,edgeOrdProperties...> worksetEdgeOrds,
726  const shards::CellTopology parentCell );
727 
767  template<typename faceTanValueType, class ...faceTanProperties,
768  typename worksetJacobianValueType, class ...worksetJacobianProperties>
769  static void
770  getPhysicalFaceTangents( Kokkos::DynRankView<faceTanValueType,faceTanProperties...> faceTanU,
771  Kokkos::DynRankView<faceTanValueType,faceTanProperties...> faceTanV,
772  const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
773  const ordinal_type worksetFaceOrd,
774  const shards::CellTopology parentCell );
775 
776 
789  template<typename faceTanValueType, class ...faceTanProperties,
790  typename worksetJacobianValueType, class ...worksetJacobianProperties,
791  typename faceOrdValueType, class ...faceOrdProperties>
792  static void
793  getPhysicalFaceTangents( Kokkos::DynRankView<faceTanValueType,faceTanProperties...> faceTanU,
794  Kokkos::DynRankView<faceTanValueType,faceTanProperties...> faceTanV,
795  const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
796  const Kokkos::DynRankView<faceOrdValueType,faceOrdProperties...> worksetFaceOrds,
797  const shards::CellTopology parentCell );
798 
799 
800 
861  template<typename sideNormalValueType, class ...sideNormalProperties,
862  typename worksetJacobianValueType, class ...worksetJacobianProperties>
863  static void
864  getPhysicalSideNormals( Kokkos::DynRankView<sideNormalValueType,sideNormalProperties...> sideNormals,
865  const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
866  const ordinal_type worksetSideOrd,
867  const shards::CellTopology parentCell );
868 
869 
880  template<typename sideNormalValueType, class ...sideNormalProperties,
881  typename worksetJacobianValueType, class ...worksetJacobianProperties,
882  typename edgeOrdValueType, class ...edgeOrdProperties>
883  static void
884  getPhysicalSideNormals( Kokkos::DynRankView<sideNormalValueType,sideNormalProperties...> sideNormals,
885  const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
886  const Kokkos::DynRankView<edgeOrdValueType,edgeOrdProperties...> worksetSideOrds,
887  const shards::CellTopology parentCell );
888 
927  template<typename faceNormalValueType, class ...faceNormalProperties,
928  typename worksetJacobianValueType, class ...worksetJacobianProperties>
929  static void
930  getPhysicalFaceNormals( Kokkos::DynRankView<faceNormalValueType,faceNormalProperties...> faceNormals,
931  const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
932  const ordinal_type worksetFaceOrd,
933  const shards::CellTopology parentCell );
934 
935 
947  template<typename faceNormalValueType, class ...faceNormalProperties,
948  typename worksetJacobianValueType, class ...worksetJacobianProperties,
949  typename faceOrdValueType, class ...faceOrdProperties>
950  static void
951  getPhysicalFaceNormals( Kokkos::DynRankView<faceNormalValueType,faceNormalProperties...> faceNormals,
952  const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
953  const Kokkos::DynRankView<faceOrdValueType,faceOrdProperties...> worksetFaceOrds,
954  const shards::CellTopology parentCell );
955 
956  //============================================================================================//
957  // //
958  // Reference-to-physical frame mapping and its inverse //
959  // //
960  //============================================================================================//
961 
998  template<typename PhysPointValueType,
999  typename RefPointValueType,
1000  typename WorksetType,
1001  typename HGradBasisPtrType>
1002  static void
1003  mapToPhysicalFrame( PhysPointValueType physPoints,
1004  const RefPointValueType refPoints,
1005  const WorksetType worksetCell,
1006  const HGradBasisPtrType basis );
1007 
1048  template<typename PhysPointViewType,
1049  typename RefPointViewType,
1050  typename WorksetCellViewType>
1051  static void
1052  mapToPhysicalFrame( PhysPointViewType physPoints,
1053  const RefPointViewType refPoints,
1054  const WorksetCellViewType worksetCell,
1055  const shards::CellTopology cellTopo ) {
1056  using nonConstRefPointValueType = typename RefPointViewType::non_const_value_type;
1057  auto basis = createHGradBasis<nonConstRefPointValueType,nonConstRefPointValueType>(cellTopo);
1058  mapToPhysicalFrame(physPoints,
1059  refPoints,
1060  worksetCell,
1061  basis);
1062  }
1063 
1064 
1115  template<typename refSubcellViewType,
1116  typename paramPointViewType>
1117  static void
1118  mapToReferenceSubcell( refSubcellViewType refSubcellPoints,
1119  const paramPointViewType paramPoints,
1120  const ordinal_type subcellDim,
1121  const ordinal_type subcellOrd,
1122  const shards::CellTopology parentCell );
1123 
1124 
1131  template<typename refSubcellViewType, typename paramPointViewType>
1132  static void
1133  mapToReferenceSubcell( refSubcellViewType refSubcellPoints,
1134  const paramPointViewType paramPoints,
1135  const typename RefSubcellParametrization<DeviceType>::ConstViewType subcellParametrization,
1136  const ordinal_type subcellOrd);
1137 
1143  template<typename refSubcellViewType, typename paramPointViewType, typename ordViewType>
1144  static void
1145  mapToReferenceSubcell( refSubcellViewType refSubcellPoints,
1146  const paramPointViewType paramPoints,
1147  const typename RefSubcellParametrization<DeviceType>::ConstViewType subcellParametrization,
1148  const ordViewType subcellOrd);
1149 
1150 
1151  //============================================================================================//
1152  // //
1153  // Physical-to-reference frame mapping and its inverse //
1154  // //
1155  //============================================================================================//
1156 
1157 
1201  template<typename refPointValueType, class ...refPointProperties,
1202  typename physPointValueType, class ...physPointProperties,
1203  typename worksetCellValueType, class ...worksetCellProperties>
1204  static void
1205  mapToReferenceFrame( Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
1206  const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1207  const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1208  const shards::CellTopology cellTopo );
1209 
1236  template<typename refPointValueType, class ...refPointProperties,
1237  typename initGuessValueType, class ...initGuessProperties,
1238  typename physPointValueType, class ...physPointProperties,
1239  typename worksetCellValueType, class ...worksetCellProperties,
1240  typename HGradBasisPtrType>
1241  static void
1242  mapToReferenceFrameInitGuess( Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
1243  const Kokkos::DynRankView<initGuessValueType,initGuessProperties...> initGuess,
1244  const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1245  const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1246  const HGradBasisPtrType basis );
1247 
1248 
1279  template<typename refPointValueType, class ...refPointProperties,
1280  typename initGuessValueType, class ...initGuessProperties,
1281  typename physPointValueType, class ...physPointProperties,
1282  typename worksetCellValueType, class ...worksetCellProperties>
1283  static void
1284  mapToReferenceFrameInitGuess( Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
1285  const Kokkos::DynRankView<initGuessValueType,initGuessProperties...> initGuess,
1286  const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1287  const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1288  const shards::CellTopology cellTopo ) {
1289  auto basis = createHGradBasis<refPointValueType,refPointValueType>(cellTopo);
1290  mapToReferenceFrameInitGuess(refPoints,
1291  initGuess,
1292  physPoints,
1293  worksetCell,
1294  basis);
1295  }
1296 
1297 
1298  //============================================================================================//
1299  // //
1300  // Control Volume Coordinates //
1301  // //
1302  //============================================================================================//
1303 
1377  template<typename subcvCoordValueType, class ...subcvCoordProperties,
1378  typename cellCoordValueType, class ...cellCoordProperties>
1379  static void
1380  getSubcvCoords( Kokkos::DynRankView<subcvCoordValueType,subcvCoordProperties...> subcvCoords,
1381  const Kokkos::DynRankView<cellCoordValueType,cellCoordProperties...> cellCoords,
1382  const shards::CellTopology primaryCell );
1383 
1384  //============================================================================================//
1385  // //
1386  // Inclusion tests //
1387  // //
1388  //============================================================================================//
1389 
1399  template<typename PointViewType>
1400  static bool
1401  checkPointInclusion( const PointViewType point,
1402  const shards::CellTopology cellTopo,
1403  const typename ScalarTraits<typename PointViewType::value_type>::scalar_type thres =
1404  threshold<typename ScalarTraits<typename PointViewType::value_type>::scalar_type>() );
1405 
1406 
1407 
1417  template<unsigned cellTopologyKey,
1418  typename OutputViewType,
1419  typename InputViewType>
1420  static void checkPointwiseInclusion( OutputViewType inCell,
1421  const InputViewType points,
1424 
1425 
1426 
1435  template<typename InCellViewType,
1436  typename PointViewType>
1437  static void checkPointwiseInclusion( InCellViewType inCell,
1438  const PointViewType points,
1439  const shards::CellTopology cellTopo,
1440  const typename ScalarTraits<typename PointViewType::value_type>::scalar_type thres =
1441  threshold<typename ScalarTraits<typename PointViewType::value_type>::scalar_type>() );
1442 
1454  template<typename inCellValueType, class ...inCellProperties,
1455  typename pointValueType, class ...pointProperties,
1456  typename cellWorksetValueType, class ...cellWorksetProperties>
1457  static void checkPointwiseInclusion( Kokkos::DynRankView<inCellValueType,inCellProperties...> inCell,
1458  const Kokkos::DynRankView<pointValueType,pointProperties...> points,
1459  const Kokkos::DynRankView<cellWorksetValueType,cellWorksetProperties...> cellWorkset,
1460  const shards::CellTopology cellTopo,
1461  const typename ScalarTraits<pointValueType>::scalar_type thres =
1462  threshold<typename ScalarTraits<pointValueType>::scalar_type>() );
1463 
1464 
1465  // //============================================================================================//
1466  // // //
1467  // // Debug //
1468  // // //
1469  // //============================================================================================//
1470 
1471 
1472  // /** \brief Prints the reference space coordinates of the vertices of the specified subcell
1473  // \param subcellDim [in] - dimension of the subcell where points are mapped to
1474  // \param subcellOrd [in] - subcell ordinal
1475  // \param parentCell [in] - cell topology of the parent cell.
1476  // */
1477  // static void printSubcellVertices(const int subcellDim,
1478  // const int subcellOrd,
1479  // const shards::CellTopology & parentCell);
1480 
1481 
1482 
1483  // /** \brief Prints the nodes of a subcell from a cell workset
1484 
1485  // */
1486  // template<class ArrayCell>
1487  // static void printWorksetSubcell(const ArrayCell & cellWorkset,
1488  // const shards::CellTopology & parentCell,
1489  // const int& pCellOrd,
1490  // const int& subcellDim,
1491  // const int& subcellOrd,
1492  // const int& fieldWidth = 3);
1493  };
1494 
1495  //============================================================================================//
1496  // //
1497  // Validation of input/output arguments for CellTools methods //
1498  // //
1499  //============================================================================================//
1500 
1507  template<typename jacobianViewType,
1508  typename PointViewType,
1509  typename worksetCellViewType>
1510  static void
1511  CellTools_setJacobianArgs( const jacobianViewType jacobian,
1512  const PointViewType points,
1513  const worksetCellViewType worksetCell,
1514  const shards::CellTopology cellTopo );
1515 
1520  template<typename jacobianInvViewType,
1521  typename jacobianViewType>
1522  static void
1523  CellTools_setJacobianInvArgs( const jacobianInvViewType jacobianInv,
1524  const jacobianViewType jacobian );
1525 
1526 
1531  template<typename jacobianDetViewType,
1532  typename jacobianViewType>
1533  static void
1534  CellTools_setJacobianDetArgs( const jacobianDetViewType jacobianDet,
1535  const jacobianViewType jacobian );
1536 
1537 
1544  template<typename physPointViewType,
1545  typename refPointViewType,
1546  typename worksetCellViewType>
1547  static void
1548  CellTools_mapToPhysicalFrameArgs( const physPointViewType physPoints,
1549  const refPointViewType refPoints,
1550  const worksetCellViewType worksetCell,
1551  const shards::CellTopology cellTopo );
1552 
1553 
1560  template<typename refPointViewType,
1561  typename physPointViewType,
1562  typename worksetCellViewType>
1563  static void
1564  CellTools_mapToReferenceFrameArgs( const refPointViewType refPoints,
1565  const physPointViewType physPoints,
1566  const worksetCellViewType worksetCell,
1567  const shards::CellTopology cellTopo );
1568 
1569 
1570 
1578  template<typename refPointViewType,
1579  typename initGuessViewType,
1580  typename physPointViewType,
1581  typename worksetCellViewType>
1582  static void
1583  CellTools_mapToReferenceFrameInitGuess( const refPointViewType refPoints,
1584  const initGuessViewType initGuess,
1585  const physPointViewType physPoints,
1586  const worksetCellViewType worksetCell,
1587  const shards::CellTopology cellTopo );
1588 
1589 }
1590 
1592 
1595 
1599 
1601 
1603 
1604 
1605 #endif
1606 
static void checkPointwiseInclusion(OutputViewType inCell, const InputViewType points, const typename ScalarTraits< typename InputViewType::value_type >::scalar_type thresh=threshold< typename ScalarTraits< typename InputViewType::value_type >::scalar_type >())
Checks every point for inclusion in the reference cell of a given topology. The points can belong to ...
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Quadrilateral cell...
Definition file for the reference to physical mappings in the Intrepid2::CellTools class...
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.
static void setJacobianDetInv(Data< PointScalar, DeviceType > &jacobianDetInv, const Data< PointScalar, DeviceType > &jacobian)
Computes reciprocals of determinants corresponding to the Jacobians in the Data container provided...
static void mapToPhysicalFrame(PhysPointViewType physPoints, const RefPointViewType refPoints, const WorksetCellViewType worksetCell, const shards::CellTopology cellTopo)
Computes F, the reference-to-physical frame map.
static void getPhysicalSideNormals(Kokkos::DynRankView< sideNormalValueType, sideNormalProperties...> sideNormals, const Kokkos::DynRankView< worksetJacobianValueType, worksetJacobianProperties...> worksetJacobians, const ordinal_type worksetSideOrd, const shards::CellTopology parentCell)
Computes non-normalized normal vectors to physical sides in a side workset .
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.
static void mapToReferenceFrameInitGuess(Kokkos::DynRankView< refPointValueType, refPointProperties...> refPoints, const Kokkos::DynRankView< initGuessValueType, initGuessProperties...> initGuess, const Kokkos::DynRankView< physPointValueType, physPointProperties...> physPoints, const Kokkos::DynRankView< worksetCellValueType, worksetCellProperties...> worksetCell, const HGradBasisPtrType basis)
Computation of , the inverse of the reference-to-physical frame map using user-supplied initial guess...
Header file for the Intrepid2::Basis_HGRAD_QUAD_C2_FEM class.
Definition file for the control volume functions of the Intrepid2::CellTools 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...
static void getReferenceCellCenter(Kokkos::DynRankView< cellCenterValueType, cellCenterProperties...> cellCenter, const shards::CellTopology cell)
Computes the Cartesian coordinates of reference cell barycenter.
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.
A stateless class for operations on cell data. Provides methods for:
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.
static void getPhysicalFaceTangents(Kokkos::DynRankView< faceTanValueType, faceTanProperties...> faceTanU, Kokkos::DynRankView< faceTanValueType, faceTanProperties...> faceTanV, const Kokkos::DynRankView< worksetJacobianValueType, worksetJacobianProperties...> worksetJacobians, const ordinal_type worksetFaceOrd, const shards::CellTopology parentCell)
Computes non-normalized tangent vector pairs to physical faces in a face workset ; (see Subcell works...
Header file for the Intrepid2::Basis_HGRAD_WEDGE_C2_FEM class.
Definition file for point inclusion functions of the Intrepid2::CellTools 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.
CellTools()=default
Default constructor.
Header function for Intrepid2::Util class and other utility functions.
Header file for the Intrepid2::Basis_HGRAD_WEDGE_C1_FEM class.
static void setJacobianDet(JacobianDetViewType jacobianDet, const JacobianViewType jacobian)
Computes the determinant of the Jacobian matrix DF of the reference-to-physical frame map F...
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.
scalar type traits
Header file for the classes: Intrepid2::RefSubcellParametrization, Intrepid2::RefCellNodes, Intrepid2::RefCellCenter.
Header file for the Intrepid2::Basis_HGRAD_HEX_C1_FEM class.
static void getReferenceEdgeTangent(RefEdgeTangentViewType refEdgeTangent, const ordinal_type edgeOrd, const shards::CellTopology parentCell)
Computes constant tangent vectors to edges of 2D or 3D reference cells.
static void mapToReferenceFrame(Kokkos::DynRankView< refPointValueType, refPointProperties...> refPoints, const Kokkos::DynRankView< physPointValueType, physPointProperties...> physPoints, const Kokkos::DynRankView< worksetCellValueType, worksetCellProperties...> worksetCell, const shards::CellTopology cellTopo)
Computes , the inverse of the reference-to-physical frame map using a default initial guess...
Definition file for the validate arguments functions of the Intrepid2::CellTools class.
Definition file for the physical to reference mappings in the Intrepid2::CellTools class...
static void getReferenceVertex(Kokkos::DynRankView< cellVertexValueType, cellVertexProperties...> cellVertex, const shards::CellTopology cell, const ordinal_type vertexOrd)
Retrieves the Cartesian coordinates of a reference cell vertex.
static bool isSupported(const unsigned cellTopoKey)
Checks if a cell topology has a reference parametrization.
static void getReferenceNode(Kokkos::DynRankView< cellNodeValueType, cellNodeProperties...> cellNode, const shards::CellTopology cell, const ordinal_type nodeOrd)
Retrieves the Cartesian coordinates of a reference cell node.
static Data< PointScalar, DeviceType > allocateJacobianInv(const Data< PointScalar, DeviceType > &jacobian)
Allocates and returns a Data container suitable for storing inverses corresponding to the Jacobians i...
static void getReferenceSideNormal(RefSideNormalViewType refSideNormal, const ordinal_type sideOrd, const shards::CellTopology parentCell)
Computes constant normal vectors to sides of 2D or 3D reference cells.
~CellTools()=default
Destructor.
static void setJacobian(JacobianViewType jacobian, const PointViewType points, const WorksetCellViewType worksetCell, const shards::CellTopology cellTopo)
Computes the Jacobian matrix DF of the reference-to-physical frame map F.
static void mapToReferenceFrameInitGuess(Kokkos::DynRankView< refPointValueType, refPointProperties...> refPoints, const Kokkos::DynRankView< initGuessValueType, initGuessProperties...> initGuess, const Kokkos::DynRankView< physPointValueType, physPointProperties...> physPoints, const Kokkos::DynRankView< worksetCellValueType, worksetCellProperties...> worksetCell, const shards::CellTopology cellTopo)
Computation of , the inverse of the reference-to-physical frame map using user-supplied initial guess...
Definition file for the Jacobian functions in the Intrepid2::CellTools class.
static void getReferenceFaceNormal(RefFaceNormalViewType refFaceNormal, const ordinal_type faceOrd, const shards::CellTopology parentCell)
Computes constant normal vectors to faces of 3D reference cell.
Definition file for node data and subcell functions of the Intrepid2::CellTools class.
Header file for the Intrepid2::Basis_HGRAD_TET_C1_FEM class.
Contains definitions of custom data types in Intrepid2.
static void setJacobian(JacobianViewType jacobian, const PointViewType points, const WorksetType worksetCell, const Teuchos::RCP< HGradBasisType > basis, const int startCell=0, const int endCell=-1)
Computes the Jacobian matrix DF of the reference-to-physical frame map F.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Quadrilateral cell...
Header file with additional documentation for the Intrepid2::CellTools class.
Header file for the Intrepid2::Basis_HGRAD_LINE_C1_FEM class.
static void setJacobianInv(JacobianInvViewType jacobianInv, const JacobianViewType jacobian)
Computes the inverse of the Jacobian matrix DF of the reference-to-physical frame map F...
static void getReferenceSubcellVertices(Kokkos::DynRankView< subcellVertexValueType, subcellVertexProperties...> subcellVertices, const ordinal_type subcellDim, const ordinal_type subcellOrd, const shards::CellTopology parentCell)
Retrieves the Cartesian coordinates of all vertices of a reference subcell.
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...
static void mapToReferenceSubcell(refSubcellViewType refSubcellPoints, const paramPointViewType paramPoints, const ordinal_type subcellDim, const ordinal_type subcellOrd, const shards::CellTopology parentCell)
Computes parameterization maps of 1- and 2-subcells of reference cells.
static void setJacobianDividedByDet(Data< PointScalar, DeviceType > &jacobianDividedByDet, const Data< PointScalar, DeviceType > &jacobian, const Data< PointScalar, DeviceType > &jacobianDetInv)
Multiplies the Jacobian with shape (C,P,D,D) by the reciprocals of the determinants, with shape (C,P), entrywise.
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.
static void getSubcvCoords(Kokkos::DynRankView< subcvCoordValueType, subcvCoordProperties...> subcvCoords, const Kokkos::DynRankView< cellCoordValueType, cellCoordProperties...> cellCoords, const shards::CellTopology primaryCell)
Computes coordinates of sub-control volumes in each primary cell.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Tetrahedron cell...
Header file for Intrepid2::RealSpaceTools class providing basic linear algebra functionality in 1D...
static void getPhysicalFaceNormals(Kokkos::DynRankView< faceNormalValueType, faceNormalProperties...> faceNormals, const Kokkos::DynRankView< worksetJacobianValueType, worksetJacobianProperties...> worksetJacobians, const ordinal_type worksetFaceOrd, const shards::CellTopology parentCell)
Computes non-normalized normal vectors to physical faces in a face workset ; (see Subcell worksets fo...
static Teuchos::RCP< Basis< DeviceType, outputValueType, pointValueType > > createHGradBasis(const shards::CellTopology cellTopo)
Generates default HGrad basis based on cell topology.
static bool hasReferenceCell(const shards::CellTopology cellTopo)
Checks if a cell topology has reference cell.
static Data< PointScalar, DeviceType > allocateJacobianDet(const Data< PointScalar, DeviceType > &jacobian)
Allocates and returns a Data container suitable for storing determinants corresponding to the Jacobia...
static void getReferenceSubcellNodes(SubcellNodeViewType subcellNodes, const ordinal_type subcellDim, const ordinal_type subcellOrd, const shards::CellTopology parentCell)
Retrieves the Cartesian coordinates of all nodes of a reference subcell.
static void mapToPhysicalFrame(PhysPointValueType physPoints, const RefPointValueType refPoints, const WorksetType worksetCell, const HGradBasisPtrType basis)
Computes F, the reference-to-physical frame map.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Hexahedron cell...
static void getPhysicalEdgeTangents(Kokkos::DynRankView< edgeTangentValueType, edgeTangentProperties...> edgeTangents, const Kokkos::DynRankView< worksetJacobianValueType, worksetJacobianProperties...> worksetJacobians, const ordinal_type worksetEdgeOrd, const shards::CellTopology parentCell)
Computes non-normalized tangent vectors to physical edges in an edge workset ; (see Subcell worksets ...
Header file for the abstract base class Intrepid2::Basis.
static void getReferenceFaceTangents(RefFaceTanViewType refFaceTanU, RefFaceTanViewType refFaceTanV, const ordinal_type faceOrd, const shards::CellTopology parentCell)
Computes pairs of constant tangent vectors to faces of a 3D reference cells.
static bool checkPointInclusion(const PointViewType point, const shards::CellTopology cellTopo, const typename ScalarTraits< typename PointViewType::value_type >::scalar_type thres=threshold< typename ScalarTraits< typename PointViewType::value_type >::scalar_type >())
Checks if a point belongs to a reference cell.