Intrepid2
Intrepid2_CellTools.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Intrepid2 Package
5 // Copyright (2007) Sandia Corporation
6 //
7 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8 // license for use of this work by or on behalf of the U.S. Government.
9 //
10 // Redistribution and use in source and binary forms, with or without
11 // modification, are permitted provided that the following conditions are
12 // met:
13 //
14 // 1. Redistributions of source code must retain the above copyright
15 // notice, this list of conditions and the following disclaimer.
16 //
17 // 2. Redistributions in binary form must reproduce the above copyright
18 // notice, this list of conditions and the following disclaimer in the
19 // documentation and/or other materials provided with the distribution.
20 //
21 // 3. Neither the name of the Corporation nor the names of the
22 // contributors may be used to endorse or promote products derived from
23 // this software without specific prior written permission.
24 //
25 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36 //
37 // Questions? Contact Kyungjoo Kim (kyukim@sandia.gov), or
38 // Mauro Perego (mperego@sandia.gov)
39 //
40 // ************************************************************************
41 // @HEADER
42 
49 #ifndef __INTREPID2_CELLTOOLS_HPP__
50 #define __INTREPID2_CELLTOOLS_HPP__
51 
52 #include "Intrepid2_ConfigDefs.hpp"
53 
54 #include "Shards_CellTopology.hpp"
55 #include "Shards_BasicTopologies.hpp"
56 
57 #include "Teuchos_RCP.hpp"
58 
59 #include "Intrepid2_Types.hpp"
60 #include "Intrepid2_Utils.hpp"
61 
63 
64 #include "Intrepid2_Basis.hpp"
65 
69 
74 
77 
81 
83 //#include "Intrepid2_HGRAD_WEDGE_I2_FEM.hpp"
84 
85 namespace Intrepid2 {
86 
87  //============================================================================================//
88  // //
89  // CellTools //
90  // //
91  //============================================================================================//
92 
103  template<typename ExecSpaceType>
104  class CellTools {
105  public:
106 
111  inline
112  static bool
113  hasReferenceCell( const shards::CellTopology cellTopo ) {
114  bool r_val = false;
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:
122 
123  case shards::Triangle<3>::key:
124  // case shards::Triangle<4>::key:
125  case shards::Triangle<6>::key:
126  // case shards::ShellTriangle<3>::key:
127  // case shards::ShellTriangle<6>::key:
128 
129  case shards::Quadrilateral<4>::key:
130  case shards::Quadrilateral<8>::key:
131  case shards::Quadrilateral<9>::key:
132  //case shards::ShellQuadrilateral<4>::key:
133  //case shards::ShellQuadrilateral<8>::key:
134  //case shards::ShellQuadrilateral<9>::key:
135 
136  case shards::Tetrahedron<4>::key:
137  // case shards::Tetrahedron<8>::key:
138  case shards::Tetrahedron<10>::key:
139  // case shards::Tetrahedron<11>::key:
140 
141  case shards::Hexahedron<8>::key:
142  case shards::Hexahedron<20>::key:
143  case shards::Hexahedron<27>::key:
144 
145  case shards::Pyramid<5>::key:
146  // case shards::Pyramid<13>::key:
147  // case shards::Pyramid<14>::key:
148 
149  case shards::Wedge<6>::key:
150  // case shards::Wedge<15>::key:
151  case shards::Wedge<18>::key:
152  r_val = true;
153  }
154  return r_val;
155  }
156 
157  private:
158 
162  template<typename outputValueType,
163  typename pointValueType>
164  static Teuchos::RCP<Basis<ExecSpaceType,outputValueType,pointValueType> >
165  createHGradBasis( const shards::CellTopology cellTopo ) {
166  Teuchos::RCP<Basis<ExecSpaceType,outputValueType,pointValueType> > r_val;
167 
168  switch (cellTopo.getKey()) {
169  case shards::Line<2>::key: r_val = Teuchos::rcp(new Basis_HGRAD_LINE_C1_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
170  case shards::Triangle<3>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TRI_C1_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
171  case shards::Quadrilateral<4>::key: r_val = Teuchos::rcp(new Basis_HGRAD_QUAD_C1_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
172  case shards::Tetrahedron<4>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TET_C1_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
173  case shards::Hexahedron<8>::key: r_val = Teuchos::rcp(new Basis_HGRAD_HEX_C1_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
174  case shards::Wedge<6>::key: r_val = Teuchos::rcp(new Basis_HGRAD_WEDGE_C1_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
175  case shards::Pyramid<5>::key: r_val = Teuchos::rcp(new Basis_HGRAD_PYR_C1_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
176 
177  case shards::Triangle<6>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TRI_C2_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
178  case shards::Quadrilateral<9>::key: r_val = Teuchos::rcp(new Basis_HGRAD_QUAD_C2_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
179  case shards::Tetrahedron<10>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TET_C2_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
180  case shards::Tetrahedron<11>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TET_COMP12_FEM<ExecSpaceType,outputValueType,pointValueType>()); break;
181  //case shards::Hexahedron<20>::key: r_val = Teuchos::rcp(new Basis_HGRAD_HEX_I2_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
182  case shards::Hexahedron<27>::key: r_val = Teuchos::rcp(new Basis_HGRAD_HEX_C2_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
183  //case shards::Wedge<15>::key: r_val = Teuchos::rcp(new Basis_HGRAD_WEDGE_I2_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
184  case shards::Wedge<18>::key: r_val = Teuchos::rcp(new Basis_HGRAD_WEDGE_C2_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
185  //case shards::Pyramid<13>::key: r_val = Teuchos::rcp(new Basis_HGRAD_PYR_I2_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
186 
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:
198  default: {
199  INTREPID2_TEST_FOR_EXCEPTION( true, std::invalid_argument,
200  ">>> ERROR (Intrepid2::CellTools::createHGradBasis): Cell topology not supported.");
201  }
202  }
203  return r_val;
204  }
205 
206  // reference nodes initialized
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];
219  };
220 
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;
233  };
234 
235  static const ReferenceNodeDataStatic refNodeDataStatic_;
236  static ReferenceNodeData refNodeData_;
237 
238  static bool isReferenceNodeDataSet_;
241  static void setReferenceNodeData();
242 
243  //============================================================================================//
244  // //
245  // Parametrization coefficients of edges and faces of reference cells //
246  // //
247  //============================================================================================//
248 
249  // static variables
253  typedef Kokkos::DynRankView<double,ExecSpaceType> subcellParamViewType;
255  subcellParamViewType dummy;
256  subcellParamViewType lineEdges; // edge maps for 2d non-standard cells; shell line and beam
257  subcellParamViewType triEdges, quadEdges; // edge maps for 2d standard cells
258  subcellParamViewType shellTriEdges, shellQuadEdges; // edge maps for 3d non-standard cells; shell tri and quad
259  subcellParamViewType tetEdges, hexEdges, pyrEdges, wedgeEdges; // edge maps for 3d standard cells
260  subcellParamViewType shellTriFaces, shellQuadFaces; // face maps for 3d non-standard cells
261  subcellParamViewType tetFaces, hexFaces, pyrFaces, wedgeFaces; // face maps for 3d standard cells
262  };
263 
264  static SubcellParamData subcellParamData_;
265 
266  static bool isSubcellParametrizationSet_;
267 
302  static void setSubcellParametrization();
303 
314  static void
315  setSubcellParametrization( subcellParamViewType &subcellParam,
316  const ordinal_type subcellDim,
317  const shards::CellTopology parentCell );
318 
319  public:
320 
323  CellTools() = default;
324 
327  ~CellTools() = default;
328 
329  //============================================================================================//
330  // //
331  // Jacobian, inverse Jacobian and Jacobian determinant //
332  // //
333  //============================================================================================//
334 
368  template<typename jacobianValueType, class ...jacobianProperties,
369  typename pointValueType, class ...pointProperties,
370  typename worksetCellValueType, class ...worksetCellProperties,
371  typename HGradBasisPtrType>
372  static void
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 );
377 
412  template<typename jacobianValueType, class ...jacobianProperties,
413  typename pointValueType, class ...pointProperties,
414  typename worksetCellValueType, class ...worksetCellProperties>
415  static void
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);
421  setJacobian(jacobian,
422  points,
423  worksetCell,
424  basis);
425  }
426 
437  template<typename jacobianInvValueType, class ...jacobianInvProperties,
438  typename jacobianValueType, class ...jacobianProperties>
439  static void
440  setJacobianInv( Kokkos::DynRankView<jacobianInvValueType,jacobianInvProperties...> jacobianInv,
441  const Kokkos::DynRankView<jacobianValueType,jacobianProperties...> jacobian );
442 
453  template<typename jacobianDetValueType, class ...jacobianDetProperties,
454  typename jacobianValueType, class ...jacobianProperties>
455  static void
456  setJacobianDet( Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
457  const Kokkos::DynRankView<jacobianValueType,jacobianProperties...> jacobian );
458 
459 
460  //============================================================================================//
461  // //
462  // Node information //
463  // //
464  //============================================================================================//
465 
466  // the node information can be used inside of kokkos functor and needs kokkos inline and
467  // exception should be an abort. for now, let's not decorate
468 
479  template<typename cellCenterValueType, class ...cellCenterProperties,
480  typename cellVertexValueType, class ...cellVertexProperties>
481  static void
482  getReferenceCellCenter( Kokkos::DynRankView<cellCenterValueType,cellCenterProperties...> cellCenter,
483  Kokkos::DynRankView<cellVertexValueType,cellVertexProperties...> cellVertex,
484  const shards::CellTopology cell );
485 
496  template<typename cellVertexValueType, class ...cellVertexProperties>
497  static void
498  getReferenceVertex( Kokkos::DynRankView<cellVertexValueType,cellVertexProperties...> cellVertex,
499  const shards::CellTopology cell,
500  const ordinal_type vertexOrd );
501 
502 
517  template<typename subcellVertexValueType, class ...subcellVertexProperties>
518  static void
519  getReferenceSubcellVertices( Kokkos::DynRankView<subcellVertexValueType,subcellVertexProperties...> subcellVertices,
520  const ordinal_type subcellDim,
521  const ordinal_type subcellOrd,
522  const shards::CellTopology parentCell );
523 
524 
525 
541  template<typename cellNodeValueType, class ...cellNodeProperties>
542  static void
543  getReferenceNode( Kokkos::DynRankView<cellNodeValueType,cellNodeProperties...> cellNode,
544  const shards::CellTopology cell,
545  const ordinal_type nodeOrd );
546 
547 
548 
562  template<typename subcellNodeValueType, class ...subcellNodeProperties>
563  static void
564  getReferenceSubcellNodes( Kokkos::DynRankView<subcellNodeValueType,subcellNodeProperties...> subcellNodes,
565  const ordinal_type subcellDim,
566  const ordinal_type subcellOrd,
567  const shards::CellTopology parentCell );
568 
594  template<typename refEdgeTangentValueType, class ...refEdgeTangentProperties>
595  static void
596  getReferenceEdgeTangent( Kokkos::DynRankView<refEdgeTangentValueType,refEdgeTangentProperties...> refEdgeTangent,
597  const ordinal_type edgeOrd,
598  const shards::CellTopology parentCell );
599 
636  template<typename refFaceTanUValueType, class ...refFaceTanUProperties,
637  typename refFaceTanVValueType, class ...refFaceTanVProperties>
638  static void
639  getReferenceFaceTangents( Kokkos::DynRankView<refFaceTanUValueType,refFaceTanUProperties...> refFaceTanU,
640  Kokkos::DynRankView<refFaceTanVValueType,refFaceTanVProperties...> refFaceTanV,
641  const ordinal_type faceOrd,
642  const shards::CellTopology parentCell );
643 
706  template<typename refSideNormalValueType, class ...refSideNormalProperties>
707  static void
708  getReferenceSideNormal( Kokkos::DynRankView<refSideNormalValueType,refSideNormalProperties...> refSideNormal,
709  const ordinal_type sideOrd,
710  const shards::CellTopology parentCell );
711 
750  template<typename refFaceNormalValueType, class ...refFaceNormalProperties>
751  static void
752  getReferenceFaceNormal( Kokkos::DynRankView<refFaceNormalValueType,refFaceNormalProperties...> refFaceNormal,
753  const ordinal_type faceOrd,
754  const shards::CellTopology parentCell );
755 
785  template<typename edgeTangentValueType, class ...edgeTangentProperties,
786  typename worksetJacobianValueType, class ...worksetJacobianProperties>
787  static void
788  getPhysicalEdgeTangents( Kokkos::DynRankView<edgeTangentValueType,edgeTangentProperties...> edgeTangents,
789  const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
790  const ordinal_type worksetEdgeOrd,
791  const shards::CellTopology parentCell );
792 
832  template<typename faceTanUValueType, class ...faceTanUProperties,
833  typename faceTanVValueType, class ...faceTanVProperties,
834  typename worksetJacobianValueType, class ...worksetJacobianProperties>
835  static void
836  getPhysicalFaceTangents( Kokkos::DynRankView<faceTanUValueType,faceTanUProperties...> faceTanU,
837  Kokkos::DynRankView<faceTanVValueType,faceTanVProperties...> faceTanV,
838  const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
839  const ordinal_type worksetFaceOrd,
840  const shards::CellTopology parentCell );
841 
902  template<typename sideNormalValueType, class ...sideNormalProperties,
903  typename worksetJacobianValueType, class ...worksetJacobianProperties>
904  static void
905  getPhysicalSideNormals( Kokkos::DynRankView<sideNormalValueType,sideNormalProperties...> sideNormals,
906  const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
907  const ordinal_type worksetSideOrd,
908  const shards::CellTopology parentCell );
909 
948  template<typename faceNormalValueType, class ...faceNormalProperties,
949  typename worksetJacobianValueType, class ...worksetJacobianProperties>
950  static void
951  getPhysicalFaceNormals( Kokkos::DynRankView<faceNormalValueType,faceNormalProperties...> faceNormals,
952  const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
953  const ordinal_type worksetFaceOrd,
954  const shards::CellTopology parentCell );
955 
956  //============================================================================================//
957  // //
958  // Reference-to-physical frame mapping and its inverse //
959  // //
960  //============================================================================================//
961 
998  template<typename physPointValueType, class ...physPointProperties,
999  typename refPointValueType, class ...refPointProperties,
1000  typename worksetCellValueType, class ...worksetCellProperties,
1001  typename HGradBasisPtrType>
1002  static void
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 );
1007 
1048  template<typename physPointValueType, class ...physPointProperties,
1049  typename refPointValueType, class ...refPointProperties,
1050  typename worksetCellValueType, class ...worksetCellProperties>
1051  static void
1052  mapToPhysicalFrame( Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
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);
1057  mapToPhysicalFrame(physPoints,
1058  refPoints,
1059  worksetCell,
1060  basis);
1061  }
1062 
1063 
1064 
1065 
1078  static void
1079  getSubcellParametrization( subcellParamViewType &subcellParam,
1080  const ordinal_type subcellDim,
1081  const shards::CellTopology parentCell );
1082 
1083 
1134  template<typename refSubcellPointValueType, class ...refSubcellPointProperties,
1135  typename paramPointValueType, class ...paramPointProperties>
1136  static void
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 );
1142 
1143 
1144  //============================================================================================//
1145  // //
1146  // Physical-to-reference frame mapping and its inverse //
1147  // //
1148  //============================================================================================//
1149 
1150 
1194  template<typename refPointValueType, class ...refPointProperties,
1195  typename physPointValueType, class ...physPointProperties,
1196  typename worksetCellValueType, class ...worksetCellProperties>
1197  static void
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 );
1202 
1229  template<typename refPointValueType, class ...refPointProperties,
1230  typename initGuessValueType, class ...initGuessProperties,
1231  typename physPointValueType, class ...physPointProperties,
1232  typename worksetCellValueType, class ...worksetCellProperties,
1233  typename HGradBasisPtrType>
1234  static void
1235  mapToReferenceFrameInitGuess( Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
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 );
1240 
1241 
1272  template<typename refPointValueType, class ...refPointProperties,
1273  typename initGuessValueType, class ...initGuessProperties,
1274  typename physPointValueType, class ...physPointProperties,
1275  typename worksetCellValueType, class ...worksetCellProperties>
1276  static void
1277  mapToReferenceFrameInitGuess( Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
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);
1283  mapToReferenceFrameInitGuess(refPoints,
1284  initGuess,
1285  physPoints,
1286  worksetCell,
1287  basis);
1288  }
1289 
1290 
1291  //============================================================================================//
1292  // //
1293  // Control Volume Coordinates //
1294  // //
1295  //============================================================================================//
1296 
1370  template<typename subcvCoordValueType, class ...subcvCoordProperties,
1371  typename cellCoordValueType, class ...cellCoordProperties>
1372  static void
1373  getSubcvCoords( Kokkos::DynRankView<subcvCoordValueType,subcvCoordProperties...> subcvCoords,
1374  const Kokkos::DynRankView<cellCoordValueType,cellCoordProperties...> cellCoords,
1375  const shards::CellTopology primaryCell );
1376 
1377  //============================================================================================//
1378  // //
1379  // Inclusion tests //
1380  // //
1381  //============================================================================================//
1382 
1392  template<typename pointValueType, class ...pointProperties>
1393  static bool
1394  checkPointInclusion( const Kokkos::DynRankView<pointValueType,pointProperties...> point,
1395  const shards::CellTopology cellTopo,
1396  const double thres = threshold() );
1397 
1398  // template<class ArrayPoint>
1399  // static int checkPointsetInclusion(const ArrayPoint & points,
1400  // const shards::CellTopology & cellTopo,
1401  // const double & threshold = INTREPID2_THRESHOLD);
1402 
1403 
1404 
1431  template<typename inCellValueType, class ...inCellProperties,
1432  typename pointValueType, class ...pointProperties>
1433  static void checkPointwiseInclusion( Kokkos::DynRankView<inCellValueType,inCellProperties...> inCell,
1434  const Kokkos::DynRankView<pointValueType,pointProperties...> points,
1435  const shards::CellTopology cellTopo,
1436  const double thres = threshold() );
1437 
1459  template<typename inCellValueType, class ...inCellProperties,
1460  typename pointValueType, class ...pointProperties,
1461  typename cellWorksetValueType, class ...cellWorksetProperties>
1462  static void checkPointwiseInclusion( Kokkos::DynRankView<inCellValueType,inCellProperties...> inCell,
1463  const Kokkos::DynRankView<pointValueType,pointProperties...> points,
1464  const Kokkos::DynRankView<cellWorksetValueType,cellWorksetProperties...> cellWorkset,
1465  const shards::CellTopology cellTopo,
1466  const double thres = threshold() );
1467 
1468 
1469  // //============================================================================================//
1470  // // //
1471  // // Debug //
1472  // // //
1473  // //============================================================================================//
1474 
1475 
1476  // /** \brief Prints the reference space coordinates of the vertices of the specified subcell
1477  // \param subcellDim [in] - dimension of the subcell where points are mapped to
1478  // \param subcellOrd [in] - subcell ordinal
1479  // \param parentCell [in] - cell topology of the parent cell.
1480  // */
1481  // static void printSubcellVertices(const int subcellDim,
1482  // const int subcellOrd,
1483  // const shards::CellTopology & parentCell);
1484 
1485 
1486 
1487  // /** \brief Prints the nodes of a subcell from a cell workset
1488 
1489  // */
1490  // template<class ArrayCell>
1491  // static void printWorksetSubcell(const ArrayCell & cellWorkset,
1492  // const shards::CellTopology & parentCell,
1493  // const int& pCellOrd,
1494  // const int& subcellDim,
1495  // const int& subcellOrd,
1496  // const int& fieldWidth = 3);
1497  };
1498 
1499  //============================================================================================//
1500  // //
1501  // Validation of input/output arguments for CellTools methods //
1502  // //
1503  //============================================================================================//
1504 
1511  template<typename jacobianViewType,
1512  typename PointViewType,
1513  typename worksetCellViewType>
1514  static void
1515  CellTools_setJacobianArgs( const jacobianViewType jacobian,
1516  const PointViewType points,
1517  const worksetCellViewType worksetCell,
1518  const shards::CellTopology cellTopo );
1519 
1524  template<typename jacobianInvViewType,
1525  typename jacobianViewType>
1526  static void
1527  CellTools_setJacobianInvArgs( const jacobianInvViewType jacobianInv,
1528  const jacobianViewType jacobian );
1529 
1530 
1535  template<typename jacobianDetViewType,
1536  typename jacobianViewType>
1537  static void
1538  CellTools_setJacobianDetArgs( const jacobianDetViewType jacobianDet,
1539  const jacobianViewType jacobian );
1540 
1541 
1548  template<typename physPointViewType,
1549  typename refPointViewType,
1550  typename worksetCellViewType>
1551  static void
1552  CellTools_mapToPhysicalFrameArgs( const physPointViewType physPoints,
1553  const refPointViewType refPoints,
1554  const worksetCellViewType worksetCell,
1555  const shards::CellTopology cellTopo );
1556 
1557 
1564  template<typename refPointViewType,
1565  typename physPointViewType,
1566  typename worksetCellViewType>
1567  static void
1568  CellTools_mapToReferenceFrameArgs( const refPointViewType refPoints,
1569  const physPointViewType physPoints,
1570  const worksetCellViewType worksetCell,
1571  const shards::CellTopology cellTopo );
1572 
1573 
1574 
1582  template<typename refPointViewType,
1583  typename initGuessViewType,
1584  typename physPointViewType,
1585  typename worksetCellViewType>
1586  static void
1587  CellTools_mapToReferenceFrameInitGuess( const refPointViewType refPoints,
1588  const initGuessViewType initGuess,
1589  const physPointViewType physPoints,
1590  const worksetCellViewType worksetCell,
1591  const shards::CellTopology cellTopo );
1592 
1593  // /** \brief Validates arguments to Intrepid2::CellTools::checkPointwiseInclusion
1594  // \param inCell [out] - rank-1 (P) array required
1595  // \param physPoints [in] - rank-2 (P,D) array required
1596  // \param cellWorkset [in] - rank-3 (C,N,D) array required
1597  // \param whichCell [in] - 0 <= whichCell < C required
1598  // \param cellTopo [in] - cell topology with a reference cell required
1599  // */
1600  // template<class ArrayIncl, class ArrayPoint, class ArrayCell>
1601  // static void
1602  // validateArguments_checkPointwiseInclusion(ArrayIncl & inCell,
1603  // const ArrayPoint & physPoints,
1604  // const ArrayCell & cellWorkset,
1605  // const int & whichCell,
1606  // const shards::CellTopology & cell);
1607 
1608 
1609 
1610 }
1611 
1613 
1616 
1619 
1622 
1624 
1625 // not yet converted ...
1627 
1628 // #include "Intrepid2_CellToolsDefDebug.hpp"
1629 
1630 
1631 #endif
1632 
Header file for the abstract base class Intrepid2::Basis.
Definition file for the control volume functions of the Intrepid2::CellTools class.
Definition file for point inclusion functions of the Intrepid2::CellTools class.
Definition file for the Jacobian functions in the Intrepid2::CellTools class.
Definition file for node data and subcell functions of the Intrepid2::CellTools class.
Definition file for the parameterization functions of the Intrepid2::CellTools class.
Definition file for the physical to reference mappings in the Intrepid2::CellTools class.
Definition file for the reference to physical mappings in the Intrepid2::CellTools class.
Definition file for the validate arguments functions of the Intrepid2::CellTools class.
Header file with additional documentation for the Intrepid2::CellTools class.
static void CellTools_mapToReferenceFrameInitGuess(const refPointViewType refPoints, const initGuessViewType initGuess, const physPointViewType physPoints, const worksetCellViewType worksetCell, const shards::CellTopology cellTopo)
Validates arguments to Intrepid2::CellTools::mapToReferenceFrame with user-defined initial guess.
Header file for the Intrepid2::Basis_HGRAD_HEX_C1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_HEX_C2_FEM class.
Header file for the Intrepid2::Basis_HGRAD_LINE_C1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_PYR_C1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_QUAD_C1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_QUAD_C2_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TET_C1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TET_C2_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TET_COMP12_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TRI_C1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TRI_C2_FEM class.
Header file for the Intrepid2::Basis_HGRAD_WEDGE_C1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_WEDGE_C2_FEM class.
Header file for Intrepid2::RealSpaceTools class providing basic linear algebra functionality in 1D,...
Contains definitions of custom data types in Intrepid2.
Header function for Intrepid2::Util class and other utility functions.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Hexahedron cell.
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 Line 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 1 on Quadrilateral cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Quadrilateral 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 2 on Tetrahedron cell.
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 1 on Triangle cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Triangle cell.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Wedge cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Wedge cell.
A stateless class for operations on cell data. Provides methods for:
static void getReferenceSideNormal(Kokkos::DynRankView< refSideNormalValueType, refSideNormalProperties... > refSideNormal, const ordinal_type sideOrd, const shards::CellTopology parentCell)
Computes constant normal vectors to sides of 2D or 3D reference cells.
static void getReferenceEdgeTangent(Kokkos::DynRankView< refEdgeTangentValueType, refEdgeTangentProperties... > refEdgeTangent, const ordinal_type edgeOrd, const shards::CellTopology parentCell)
Computes constant tangent vectors to edges of 2D or 3D reference cells.
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 void setJacobianDet(Kokkos::DynRankView< jacobianDetValueType, jacobianDetProperties... > jacobianDet, const Kokkos::DynRankView< jacobianValueType, jacobianProperties... > jacobian)
Computes the determinant of the Jacobian matrix DF of the reference-to-physical frame map F.
static void mapToReferenceSubcell(Kokkos::DynRankView< refSubcellPointValueType, refSubcellPointProperties... > refSubcellPoints, const Kokkos::DynRankView< paramPointValueType, paramPointProperties... > 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 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...
CellTools()=default
Default constructor.
static void mapToPhysicalFrame(Kokkos::DynRankView< physPointValueType, physPointProperties... > physPoints, const Kokkos::DynRankView< refPointValueType, refPointProperties... > refPoints, const Kokkos::DynRankView< worksetCellValueType, worksetCellProperties... > worksetCell, const HGradBasisPtrType basis)
Computes F, the reference-to-physical frame map.
static bool checkPointInclusion(const Kokkos::DynRankView< pointValueType, pointProperties... > point, const shards::CellTopology cellTopo, const double thres=threshold())
Checks if a point belongs to a reference cell.
static void getReferenceSubcellNodes(Kokkos::DynRankView< subcellNodeValueType, subcellNodeProperties... > 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(Kokkos::DynRankView< physPointValueType, physPointProperties... > physPoints, const Kokkos::DynRankView< refPointValueType, refPointProperties... > refPoints, const Kokkos::DynRankView< worksetCellValueType, worksetCellProperties... > worksetCell, const shards::CellTopology cellTopo)
Computes F, the reference-to-physical frame map.
static void setReferenceNodeData()
Set reference node coordinates for supported topologies.
static void checkPointwiseInclusion(Kokkos::DynRankView< inCellValueType, inCellProperties... > inCell, const Kokkos::DynRankView< pointValueType, pointProperties... > points, const shards::CellTopology cellTopo, const double thres=threshold())
Checks every point in a set for inclusion in a reference cell.
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.
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.
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 .
static void setJacobian(Kokkos::DynRankView< jacobianValueType, jacobianProperties... > jacobian, const Kokkos::DynRankView< pointValueType, pointProperties... > points, const Kokkos::DynRankView< worksetCellValueType, worksetCellProperties... > worksetCell, const HGradBasisPtrType basis)
Computes the Jacobian matrix DF of the reference-to-physical frame map F.
static bool hasReferenceCell(const shards::CellTopology cellTopo)
Checks if a cell topology has reference cell.
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 void setSubcellParametrization()
Defines orientation-preserving parametrizations of reference edges and faces of cell topologies with ...
static void getPhysicalFaceTangents(Kokkos::DynRankView< faceTanUValueType, faceTanUProperties... > faceTanU, Kokkos::DynRankView< faceTanVValueType, faceTanVProperties... > 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...
static void getSubcellParametrization(subcellParamViewType &subcellParam, const ordinal_type subcellDim, const shards::CellTopology parentCell)
Returns array with the coefficients of the parametrization maps for the edges or faces of a reference...
static void getReferenceFaceTangents(Kokkos::DynRankView< refFaceTanUValueType, refFaceTanUProperties... > refFaceTanU, Kokkos::DynRankView< refFaceTanVValueType, refFaceTanVProperties... > refFaceTanV, const ordinal_type faceOrd, const shards::CellTopology parentCell)
Computes pairs of constant tangent vectors to faces of a 3D reference cells.
static void setJacobianInv(Kokkos::DynRankView< jacobianInvValueType, jacobianInvProperties... > jacobianInv, const Kokkos::DynRankView< jacobianValueType, jacobianProperties... > jacobian)
Computes the inverse of the Jacobian matrix DF of the reference-to-physical frame map F.
~CellTools()=default
Destructor.
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 ...
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...
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.
static Teuchos::RCP< Basis< ExecSpaceType, outputValueType, pointValueType > > createHGradBasis(const shards::CellTopology cellTopo)
Generates default HGrad basis based on cell topology.
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...
static void getReferenceCellCenter(Kokkos::DynRankView< cellCenterValueType, cellCenterProperties... > cellCenter, Kokkos::DynRankView< cellVertexValueType, cellVertexProperties... > cellVertex, const shards::CellTopology cell)
Computes the Cartesian coordinates of reference cell center.
static void setJacobian(Kokkos::DynRankView< jacobianValueType, jacobianProperties... > jacobian, const Kokkos::DynRankView< pointValueType, pointProperties... > points, const Kokkos::DynRankView< worksetCellValueType, worksetCellProperties... > worksetCell, const shards::CellTopology cellTopo)
Computes the Jacobian matrix DF of the reference-to-physical frame map F.
static void getReferenceFaceNormal(Kokkos::DynRankView< refFaceNormalValueType, refFaceNormalProperties... > refFaceNormal, const ordinal_type faceOrd, const shards::CellTopology parentCell)
Computes constant normal vectors to faces of 3D reference cell.
Reference node containers for each supported topology.
Reference node data for each supported topology.
Parametrization coefficients of edges and faces of reference cells.