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  typedef Kokkos::DynRankView<double,ExecSpaceType> subcellParamViewType;
158 
159  private:
160 
164  template<typename outputValueType,
165  typename pointValueType>
166  static Teuchos::RCP<Basis<ExecSpaceType,outputValueType,pointValueType> >
167  createHGradBasis( const shards::CellTopology cellTopo ) {
168  Teuchos::RCP<Basis<ExecSpaceType,outputValueType,pointValueType> > r_val;
169 
170  switch (cellTopo.getKey()) {
171  case shards::Line<2>::key: r_val = Teuchos::rcp(new Basis_HGRAD_LINE_C1_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
172  case shards::Triangle<3>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TRI_C1_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
173  case shards::Quadrilateral<4>::key: r_val = Teuchos::rcp(new Basis_HGRAD_QUAD_C1_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
174  case shards::Tetrahedron<4>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TET_C1_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
175  case shards::Hexahedron<8>::key: r_val = Teuchos::rcp(new Basis_HGRAD_HEX_C1_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
176  case shards::Wedge<6>::key: r_val = Teuchos::rcp(new Basis_HGRAD_WEDGE_C1_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
177  case shards::Pyramid<5>::key: r_val = Teuchos::rcp(new Basis_HGRAD_PYR_C1_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
178 
179  case shards::Triangle<6>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TRI_C2_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
180  case shards::Quadrilateral<9>::key: r_val = Teuchos::rcp(new Basis_HGRAD_QUAD_C2_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
181  case shards::Tetrahedron<10>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TET_C2_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
182  case shards::Tetrahedron<11>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TET_COMP12_FEM<ExecSpaceType,outputValueType,pointValueType>()); break;
183  //case shards::Hexahedron<20>::key: r_val = Teuchos::rcp(new Basis_HGRAD_HEX_I2_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
184  case shards::Hexahedron<27>::key: r_val = Teuchos::rcp(new Basis_HGRAD_HEX_C2_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
185  //case shards::Wedge<15>::key: r_val = Teuchos::rcp(new Basis_HGRAD_WEDGE_I2_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
186  case shards::Wedge<18>::key: r_val = Teuchos::rcp(new Basis_HGRAD_WEDGE_C2_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
187  //case shards::Pyramid<13>::key: r_val = Teuchos::rcp(new Basis_HGRAD_PYR_I2_FEM <ExecSpaceType,outputValueType,pointValueType>()); break;
188 
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:
200  default: {
201  INTREPID2_TEST_FOR_EXCEPTION( true, std::invalid_argument,
202  ">>> ERROR (Intrepid2::CellTools::createHGradBasis): Cell topology not supported.");
203  }
204  }
205  return r_val;
206  }
207 
208  // reference nodes initialized
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];
221  };
222 
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;
235  };
236 
237  static const ReferenceNodeDataStatic refNodeDataStatic_;
238  static ReferenceNodeData refNodeData_;
239 
240  static bool isReferenceNodeDataSet_;
243  static void setReferenceNodeData();
244 
245  //============================================================================================//
246  // //
247  // Parametrization coefficients of edges and faces of reference cells //
248  // //
249  //============================================================================================//
250 
251  // static variables
256  subcellParamViewType dummy;
257  subcellParamViewType lineEdges; // edge maps for 2d non-standard cells; shell line and beam
258  subcellParamViewType triEdges, quadEdges; // edge maps for 2d standard cells
259  subcellParamViewType shellTriEdges, shellQuadEdges; // edge maps for 3d non-standard cells; shell tri and quad
260  subcellParamViewType tetEdges, hexEdges, pyrEdges, wedgeEdges; // edge maps for 3d standard cells
261  subcellParamViewType shellTriFaces, shellQuadFaces; // face maps for 3d non-standard cells
262  subcellParamViewType tetFaces, hexFaces, pyrFaces, wedgeFaces; // face maps for 3d standard cells
263  };
264 
265  static SubcellParamData subcellParamData_;
266 
267  static bool isSubcellParametrizationSet_;
268 
269 
280  static void
281  setSubcellParametrization( subcellParamViewType &subcellParam,
282  const ordinal_type subcellDim,
283  const shards::CellTopology parentCell );
284 
285  public:
286 
321  static void setSubcellParametrization();
322 
325  CellTools() = default;
326 
329  ~CellTools() = default;
330 
331  //============================================================================================//
332  // //
333  // Jacobian, inverse Jacobian and Jacobian determinant //
334  // //
335  //============================================================================================//
336 
370  template<typename jacobianValueType, class ...jacobianProperties,
371  typename pointValueType, class ...pointProperties,
372  typename worksetCellValueType, class ...worksetCellProperties,
373  typename HGradBasisPtrType>
374  static void
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 );
379 
414  template<typename jacobianValueType, class ...jacobianProperties,
415  typename pointValueType, class ...pointProperties,
416  typename worksetCellValueType, class ...worksetCellProperties>
417  static void
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);
423  setJacobian(jacobian,
424  points,
425  worksetCell,
426  basis);
427  }
428 
439  template<typename jacobianInvValueType, class ...jacobianInvProperties,
440  typename jacobianValueType, class ...jacobianProperties>
441  static void
442  setJacobianInv( Kokkos::DynRankView<jacobianInvValueType,jacobianInvProperties...> jacobianInv,
443  const Kokkos::DynRankView<jacobianValueType,jacobianProperties...> jacobian );
444 
455  template<typename jacobianDetValueType, class ...jacobianDetProperties,
456  typename jacobianValueType, class ...jacobianProperties>
457  static void
458  setJacobianDet( Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
459  const Kokkos::DynRankView<jacobianValueType,jacobianProperties...> jacobian );
460 
461 
462  //============================================================================================//
463  // //
464  // Node information //
465  // //
466  //============================================================================================//
467 
468  // the node information can be used inside of kokkos functor and needs kokkos inline and
469  // exception should be an abort. for now, let's not decorate
470 
481  template<typename cellCenterValueType, class ...cellCenterProperties,
482  typename cellVertexValueType, class ...cellVertexProperties>
483  static void
484  getReferenceCellCenter( Kokkos::DynRankView<cellCenterValueType,cellCenterProperties...> cellCenter,
485  Kokkos::DynRankView<cellVertexValueType,cellVertexProperties...> cellVertex,
486  const shards::CellTopology cell );
487 
498  template<typename cellVertexValueType, class ...cellVertexProperties>
499  static void
500  getReferenceVertex( Kokkos::DynRankView<cellVertexValueType,cellVertexProperties...> cellVertex,
501  const shards::CellTopology cell,
502  const ordinal_type vertexOrd );
503 
504 
519  template<typename subcellVertexValueType, class ...subcellVertexProperties>
520  static void
521  getReferenceSubcellVertices( Kokkos::DynRankView<subcellVertexValueType,subcellVertexProperties...> subcellVertices,
522  const ordinal_type subcellDim,
523  const ordinal_type subcellOrd,
524  const shards::CellTopology parentCell );
525 
526 
527 
543  template<typename cellNodeValueType, class ...cellNodeProperties>
544  static void
545  getReferenceNode( Kokkos::DynRankView<cellNodeValueType,cellNodeProperties...> cellNode,
546  const shards::CellTopology cell,
547  const ordinal_type nodeOrd );
548 
549 
550 
564  template<typename subcellNodeValueType, class ...subcellNodeProperties>
565  static void
566  getReferenceSubcellNodes( Kokkos::DynRankView<subcellNodeValueType,subcellNodeProperties...> subcellNodes,
567  const ordinal_type subcellDim,
568  const ordinal_type subcellOrd,
569  const shards::CellTopology parentCell );
570 
596  template<typename refEdgeTangentValueType, class ...refEdgeTangentProperties>
597  static void
598  getReferenceEdgeTangent( Kokkos::DynRankView<refEdgeTangentValueType,refEdgeTangentProperties...> refEdgeTangent,
599  const ordinal_type edgeOrd,
600  const shards::CellTopology parentCell );
601 
638  template<typename refFaceTanUValueType, class ...refFaceTanUProperties,
639  typename refFaceTanVValueType, class ...refFaceTanVProperties>
640  static void
641  getReferenceFaceTangents( Kokkos::DynRankView<refFaceTanUValueType,refFaceTanUProperties...> refFaceTanU,
642  Kokkos::DynRankView<refFaceTanVValueType,refFaceTanVProperties...> refFaceTanV,
643  const ordinal_type faceOrd,
644  const shards::CellTopology parentCell );
645 
708  template<typename refSideNormalValueType, class ...refSideNormalProperties>
709  static void
710  getReferenceSideNormal( Kokkos::DynRankView<refSideNormalValueType,refSideNormalProperties...> refSideNormal,
711  const ordinal_type sideOrd,
712  const shards::CellTopology parentCell );
713 
752  template<typename refFaceNormalValueType, class ...refFaceNormalProperties>
753  static void
754  getReferenceFaceNormal( Kokkos::DynRankView<refFaceNormalValueType,refFaceNormalProperties...> refFaceNormal,
755  const ordinal_type faceOrd,
756  const shards::CellTopology parentCell );
757 
787  template<typename edgeTangentValueType, class ...edgeTangentProperties,
788  typename worksetJacobianValueType, class ...worksetJacobianProperties>
789  static void
790  getPhysicalEdgeTangents( Kokkos::DynRankView<edgeTangentValueType,edgeTangentProperties...> edgeTangents,
791  const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
792  const ordinal_type worksetEdgeOrd,
793  const shards::CellTopology parentCell );
794 
834  template<typename faceTanUValueType, class ...faceTanUProperties,
835  typename faceTanVValueType, class ...faceTanVProperties,
836  typename worksetJacobianValueType, class ...worksetJacobianProperties>
837  static void
838  getPhysicalFaceTangents( Kokkos::DynRankView<faceTanUValueType,faceTanUProperties...> faceTanU,
839  Kokkos::DynRankView<faceTanVValueType,faceTanVProperties...> faceTanV,
840  const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
841  const ordinal_type worksetFaceOrd,
842  const shards::CellTopology parentCell );
843 
904  template<typename sideNormalValueType, class ...sideNormalProperties,
905  typename worksetJacobianValueType, class ...worksetJacobianProperties>
906  static void
907  getPhysicalSideNormals( Kokkos::DynRankView<sideNormalValueType,sideNormalProperties...> sideNormals,
908  const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
909  const ordinal_type worksetSideOrd,
910  const shards::CellTopology parentCell );
911 
950  template<typename faceNormalValueType, class ...faceNormalProperties,
951  typename worksetJacobianValueType, class ...worksetJacobianProperties>
952  static void
953  getPhysicalFaceNormals( Kokkos::DynRankView<faceNormalValueType,faceNormalProperties...> faceNormals,
954  const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
955  const ordinal_type worksetFaceOrd,
956  const shards::CellTopology parentCell );
957 
958  //============================================================================================//
959  // //
960  // Reference-to-physical frame mapping and its inverse //
961  // //
962  //============================================================================================//
963 
1000  template<typename physPointValueType, class ...physPointProperties,
1001  typename refPointValueType, class ...refPointProperties,
1002  typename worksetCellValueType, class ...worksetCellProperties,
1003  typename HGradBasisPtrType>
1004  static void
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 );
1009 
1050  template<typename physPointValueType, class ...physPointProperties,
1051  typename refPointValueType, class ...refPointProperties,
1052  typename worksetCellValueType, class ...worksetCellProperties>
1053  static void
1054  mapToPhysicalFrame( Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
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);
1059  mapToPhysicalFrame(physPoints,
1060  refPoints,
1061  worksetCell,
1062  basis);
1063  }
1064 
1065 
1066 
1067 
1080  static void
1081  getSubcellParametrization( subcellParamViewType &subcellParam,
1082  const ordinal_type subcellDim,
1083  const shards::CellTopology parentCell );
1084 
1085 
1136  template<typename refSubcellPointValueType, class ...refSubcellPointProperties,
1137  typename paramPointValueType, class ...paramPointProperties>
1138  static void
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 );
1144 
1145 
1156  template<typename refSubcellPointValueType, class ...refSubcellPointProperties,
1157  typename paramPointValueType, class ...paramPointProperties>
1158  static void
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);
1166 
1167 
1168  //============================================================================================//
1169  // //
1170  // Physical-to-reference frame mapping and its inverse //
1171  // //
1172  //============================================================================================//
1173 
1174 
1218  template<typename refPointValueType, class ...refPointProperties,
1219  typename physPointValueType, class ...physPointProperties,
1220  typename worksetCellValueType, class ...worksetCellProperties>
1221  static void
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 );
1226 
1253  template<typename refPointValueType, class ...refPointProperties,
1254  typename initGuessValueType, class ...initGuessProperties,
1255  typename physPointValueType, class ...physPointProperties,
1256  typename worksetCellValueType, class ...worksetCellProperties,
1257  typename HGradBasisPtrType>
1258  static void
1259  mapToReferenceFrameInitGuess( Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
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 );
1264 
1265 
1296  template<typename refPointValueType, class ...refPointProperties,
1297  typename initGuessValueType, class ...initGuessProperties,
1298  typename physPointValueType, class ...physPointProperties,
1299  typename worksetCellValueType, class ...worksetCellProperties>
1300  static void
1301  mapToReferenceFrameInitGuess( Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
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);
1307  mapToReferenceFrameInitGuess(refPoints,
1308  initGuess,
1309  physPoints,
1310  worksetCell,
1311  basis);
1312  }
1313 
1314 
1315  //============================================================================================//
1316  // //
1317  // Control Volume Coordinates //
1318  // //
1319  //============================================================================================//
1320 
1394  template<typename subcvCoordValueType, class ...subcvCoordProperties,
1395  typename cellCoordValueType, class ...cellCoordProperties>
1396  static void
1397  getSubcvCoords( Kokkos::DynRankView<subcvCoordValueType,subcvCoordProperties...> subcvCoords,
1398  const Kokkos::DynRankView<cellCoordValueType,cellCoordProperties...> cellCoords,
1399  const shards::CellTopology primaryCell );
1400 
1401  //============================================================================================//
1402  // //
1403  // Inclusion tests //
1404  // //
1405  //============================================================================================//
1406 
1416  template<typename pointValueType, class ...pointProperties>
1417  static bool
1418  checkPointInclusion( const Kokkos::DynRankView<pointValueType,pointProperties...> point,
1419  const shards::CellTopology cellTopo,
1420  const double thres = threshold() );
1421 
1422  // template<class ArrayPoint>
1423  // static int checkPointsetInclusion(const ArrayPoint & points,
1424  // const shards::CellTopology & cellTopo,
1425  // const double & threshold = INTREPID2_THRESHOLD);
1426 
1427 
1428 
1455  template<typename inCellValueType, class ...inCellProperties,
1456  typename pointValueType, class ...pointProperties>
1457  static void checkPointwiseInclusion( Kokkos::DynRankView<inCellValueType,inCellProperties...> inCell,
1458  const Kokkos::DynRankView<pointValueType,pointProperties...> points,
1459  const shards::CellTopology cellTopo,
1460  const double thres = threshold() );
1461 
1483  template<typename inCellValueType, class ...inCellProperties,
1484  typename pointValueType, class ...pointProperties,
1485  typename cellWorksetValueType, class ...cellWorksetProperties>
1486  static void checkPointwiseInclusion( Kokkos::DynRankView<inCellValueType,inCellProperties...> inCell,
1487  const Kokkos::DynRankView<pointValueType,pointProperties...> points,
1488  const Kokkos::DynRankView<cellWorksetValueType,cellWorksetProperties...> cellWorkset,
1489  const shards::CellTopology cellTopo,
1490  const double thres = threshold() );
1491 
1492 
1493  // //============================================================================================//
1494  // // //
1495  // // Debug //
1496  // // //
1497  // //============================================================================================//
1498 
1499 
1500  // /** \brief Prints the reference space coordinates of the vertices of the specified subcell
1501  // \param subcellDim [in] - dimension of the subcell where points are mapped to
1502  // \param subcellOrd [in] - subcell ordinal
1503  // \param parentCell [in] - cell topology of the parent cell.
1504  // */
1505  // static void printSubcellVertices(const int subcellDim,
1506  // const int subcellOrd,
1507  // const shards::CellTopology & parentCell);
1508 
1509 
1510 
1511  // /** \brief Prints the nodes of a subcell from a cell workset
1512 
1513  // */
1514  // template<class ArrayCell>
1515  // static void printWorksetSubcell(const ArrayCell & cellWorkset,
1516  // const shards::CellTopology & parentCell,
1517  // const int& pCellOrd,
1518  // const int& subcellDim,
1519  // const int& subcellOrd,
1520  // const int& fieldWidth = 3);
1521  };
1522 
1523  //============================================================================================//
1524  // //
1525  // Validation of input/output arguments for CellTools methods //
1526  // //
1527  //============================================================================================//
1528 
1535  template<typename jacobianViewType,
1536  typename PointViewType,
1537  typename worksetCellViewType>
1538  static void
1539  CellTools_setJacobianArgs( const jacobianViewType jacobian,
1540  const PointViewType points,
1541  const worksetCellViewType worksetCell,
1542  const shards::CellTopology cellTopo );
1543 
1548  template<typename jacobianInvViewType,
1549  typename jacobianViewType>
1550  static void
1551  CellTools_setJacobianInvArgs( const jacobianInvViewType jacobianInv,
1552  const jacobianViewType jacobian );
1553 
1554 
1559  template<typename jacobianDetViewType,
1560  typename jacobianViewType>
1561  static void
1562  CellTools_setJacobianDetArgs( const jacobianDetViewType jacobianDet,
1563  const jacobianViewType jacobian );
1564 
1565 
1572  template<typename physPointViewType,
1573  typename refPointViewType,
1574  typename worksetCellViewType>
1575  static void
1576  CellTools_mapToPhysicalFrameArgs( const physPointViewType physPoints,
1577  const refPointViewType refPoints,
1578  const worksetCellViewType worksetCell,
1579  const shards::CellTopology cellTopo );
1580 
1581 
1588  template<typename refPointViewType,
1589  typename physPointViewType,
1590  typename worksetCellViewType>
1591  static void
1592  CellTools_mapToReferenceFrameArgs( const refPointViewType refPoints,
1593  const physPointViewType physPoints,
1594  const worksetCellViewType worksetCell,
1595  const shards::CellTopology cellTopo );
1596 
1597 
1598 
1606  template<typename refPointViewType,
1607  typename initGuessViewType,
1608  typename physPointViewType,
1609  typename worksetCellViewType>
1610  static void
1611  CellTools_mapToReferenceFrameInitGuess( const refPointViewType refPoints,
1612  const initGuessViewType initGuess,
1613  const physPointViewType physPoints,
1614  const worksetCellViewType worksetCell,
1615  const shards::CellTopology cellTopo );
1616 
1617  // /** \brief Validates arguments to Intrepid2::CellTools::checkPointwiseInclusion
1618  // \param inCell [out] - rank-1 (P) array required
1619  // \param physPoints [in] - rank-2 (P,D) array required
1620  // \param cellWorkset [in] - rank-3 (C,N,D) array required
1621  // \param whichCell [in] - 0 <= whichCell < C required
1622  // \param cellTopo [in] - cell topology with a reference cell required
1623  // */
1624  // template<class ArrayIncl, class ArrayPoint, class ArrayCell>
1625  // static void
1626  // validateArguments_checkPointwiseInclusion(ArrayIncl & inCell,
1627  // const ArrayPoint & physPoints,
1628  // const ArrayCell & cellWorkset,
1629  // const int & whichCell,
1630  // const shards::CellTopology & cell);
1631 
1632 
1633 
1634 }
1635 
1637 
1640 
1643 
1646 
1648 
1649 // not yet converted ...
1651 
1652 // #include "Intrepid2_CellToolsDefDebug.hpp"
1653 
1654 
1655 #endif
1656 
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.