Intrepid
Intrepid_CellTools.hpp
Go to the documentation of this file.
1 #ifndef INTREPID_CELTOOLS_HPP
2 #define INTREPID_CELTOOLS_HPP
3 
4 // @HEADER
5 // ************************************************************************
6 //
7 // Intrepid Package
8 // Copyright (2007) Sandia Corporation
9 //
10 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
11 // license for use of this work by or on behalf of the U.S. Government.
12 //
13 // Redistribution and use in source and binary forms, with or without
14 // modification, are permitted provided that the following conditions are
15 // met:
16 //
17 // 1. Redistributions of source code must retain the above copyright
18 // notice, this list of conditions and the following disclaimer.
19 //
20 // 2. Redistributions in binary form must reproduce the above copyright
21 // notice, this list of conditions and the following disclaimer in the
22 // documentation and/or other materials provided with the distribution.
23 //
24 // 3. Neither the name of the Corporation nor the names of the
25 // contributors may be used to endorse or promote products derived from
26 // this software without specific prior written permission.
27 //
28 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
29 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
30 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
31 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
32 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
33 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
34 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
35 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
36 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
37 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
38 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
39 //
40 // Questions? Contact Pavel Bochev (pbboche@sandia.gov)
41 // Denis Ridzal (dridzal@sandia.gov), or
42 // Kara Peterson (kjpeter@sandia.gov)
43 //
44 // ************************************************************************
45 // @HEADER
46 
52 #ifndef INTREPID_CELLTOOLS_HPP
53 #define INTREPID_CELLTOOLS_HPP
54 
55 
58 #include "Intrepid_ConfigDefs.hpp"
59 #include "Intrepid_Types.hpp"
60 #include "Intrepid_Utils.hpp"
61 #include "Intrepid_Basis.hpp"
63 #include "Intrepid_HGRAD_QUAD_C1_FEM.hpp"
64 #include "Intrepid_HGRAD_TET_C1_FEM.hpp"
65 #include "Intrepid_HGRAD_WEDGE_C1_FEM.hpp"
66 #include "Intrepid_HGRAD_PYR_C1_FEM.hpp"
67 #include "Intrepid_HGRAD_HEX_C1_FEM.hpp"
68 
70 
72 #include "Intrepid_HGRAD_QUAD_C2_FEM.hpp"
80 
81 #include "Shards_CellTopology.hpp"
82 #include "Shards_BasicTopologies.hpp"
83 
84 #include "Teuchos_Assert.hpp"
85 #include "Teuchos_RCP.hpp"
86 
87 #include <Intrepid_KokkosRank.hpp>
88 #ifdef INTREPID_OLD_KOKKOS_CODE
89 #include "Kokkos_Core.hpp"
90 #endif
91 namespace Intrepid {
92 
93 //nn
94  //============================================================================================//
95  // //
96  // CellTools //
97  // //
98  //============================================================================================//
99 
110 template<class Scalar>
111 class CellTools {
112 private:
113 
114  //============================================================================================//
115  // //
116  // Parametrization coefficients of edges and faces of reference cells //
117  // //
118  //============================================================================================//
119 
120 
133  static const FieldContainer<double>& getSubcellParametrization(const int subcellDim,
134  const shards::CellTopology& parentCell);
135 
136 
137 
177  static void setSubcellParametrization(FieldContainer<double>& subcellParametrization,
178  const int subcellDim,
179  const shards::CellTopology& parentCell);
180 
181  //============================================================================================//
182  // //
183  // Validation of input/output arguments for CellTools methods //
184  // //
185  //============================================================================================//
186 
194  template<class ArrayJac, class ArrayPoint, class ArrayCell>
195  static void validateArguments_setJacobian(const ArrayJac & jacobian,
196  const ArrayPoint & points,
197  const ArrayCell & cellWorkset,
198  const int & whichCell,
199  const shards::CellTopology & cellTopo);
200 
201 
202 
207  template<class ArrayJacInv, class ArrayJac>
208  static void validateArguments_setJacobianInv(const ArrayJacInv & jacobianInv,
209  const ArrayJac & jacobian);
210 
211 
212 
217  template<class ArrayJacDet, class ArrayJac>
218  static void validateArguments_setJacobianDetArgs(const ArrayJacDet & jacobianDet,
219  const ArrayJac & jacobian);
220 
221 
222 
230  template<class ArrayPhysPoint, class ArrayRefPoint, class ArrayCell>
231  static void validateArguments_mapToPhysicalFrame(const ArrayPhysPoint & physPoints,
232  const ArrayRefPoint & refPoints,
233  const ArrayCell & cellWorkset,
234  const shards::CellTopology & cellTopo,
235  const int& whichCell);
236 
237 
238 
246  template<class ArrayRefPoint, class ArrayPhysPoint, class ArrayCell>
247  static void validateArguments_mapToReferenceFrame(const ArrayRefPoint & refPoints,
248  const ArrayPhysPoint & physPoints,
249  const ArrayCell & cellWorkset,
250  const shards::CellTopology & cellTopo,
251  const int& whichCell);
252 
253 
254 
263  template<class ArrayRefPoint, class ArrayInitGuess, class ArrayPhysPoint, class ArrayCell>
264  static void validateArguments_mapToReferenceFrame(const ArrayRefPoint & refPoints,
265  const ArrayInitGuess & initGuess,
266  const ArrayPhysPoint & physPoints,
267  const ArrayCell & cellWorkset,
268  const shards::CellTopology & cellTopo,
269  const int& whichCell);
270 
271 
272 
280  template<class ArrayIncl, class ArrayPoint, class ArrayCell>
281  static void validateArguments_checkPointwiseInclusion(ArrayIncl & inCell,
282  const ArrayPoint & physPoints,
283  const ArrayCell & cellWorkset,
284  const int & whichCell,
285  const shards::CellTopology & cell);
286 public:
287 
290  CellTools(){ };
291 
292 
296 
297  //============================================================================================//
298  // //
299  // Jacobian, inverse Jacobian and Jacobian determinant //
300  // //
301  //============================================================================================//
302 
350 /*
351  #ifdef INTREPID_OLD_KOKKOS_CODE
352 
353 
354 template<class ArrayJac, class ArrayPoint, class ArrayCell, bool typecheck>
355  struct setJacobianTempSpecKokkos;
356 
357 
358  template<class Scalar1,class Scalar2,class Scalar3,class Layout,class MemorySpace>
359  static void setJacobianTemp(Kokkos::View<Scalar1,Layout,MemorySpace> & jacobian,
360  const Kokkos::View<Scalar2,Layout,MemorySpace> & points,
361  const Kokkos::View<Scalar3,Layout,MemorySpace> & cellWorkset,
362  const shards::CellTopology & cellTopo,
363  const int & whichCell = -1);
364 
365  #endif */
366  template<class ArrayJac, class ArrayPoint, class ArrayCell, bool typecheck>
368 
369 
370  /*template<class ArrayJac, class ArrayPoint, class ArrayCell>
371  static void setJacobianTemp(ArrayJac & jacobian,
372  const ArrayPoint & points,
373  const ArrayCell & cellWorkset,
374  const shards::CellTopology & cellTopo,
375  const int & whichCell = -1);*/
376  template<class ArrayJac, class ArrayPoint, class ArrayCell>
377  static void setJacobian(ArrayJac & jacobian,
378  const ArrayPoint & points,
379  const ArrayCell & cellWorkset,
380  const shards::CellTopology & cellTopo,
381  const int & whichCell = -1);
382 
383  template<class ArrayJac, class ArrayPoint, class ArrayCell>
384  static void setJacobian(ArrayJac & jacobian,
385  const ArrayPoint & points,
386  const ArrayCell & cellWorkset,
387  const Teuchos::RCP< Basis< Scalar, FieldContainer<Scalar> > > HGRAD_Basis,
388  const int & whichCell = -1);
389 
390 
403  template<class ArrayJacInv, class ArrayJac>
404  static void setJacobianInv(ArrayJacInv & jacobianInv,
405  const ArrayJac & jacobian);
406 
407  /* template<class ArrayJacInv, class ArrayJac>
408  static void setJacobianInvTemp(ArrayJacInv & jacobianInv,
409  const ArrayJac & jacobian);
410  */
411 
424  template<class ArrayJacDet, class ArrayJac>
425  static void setJacobianDet(ArrayJacDet & jacobianDet,
426  const ArrayJac & jacobian);
427 
428  /* template<class ArrayJacDet, class ArrayJac>
429  static void setJacobianDetTemp(ArrayJacDet & jacobianDet,
430  const ArrayJac & jacobian);*/
431 
432  //============================================================================================//
433  // //
434  // Reference-to-physical frame mapping and its inverse //
435  // //
436  //============================================================================================//
437 
493  template<class ArrayPhysPoint, class ArrayRefPoint, class ArrayCell>
494  static void mapToPhysicalFrame(ArrayPhysPoint & physPoints,
495  const ArrayRefPoint & refPoints,
496  const ArrayCell & cellWorkset,
497  const shards::CellTopology & cellTopo,
498  const int & whichCell = -1);
499 
500  template<class ArrayPhysPoint, class ArrayRefPoint, class ArrayCell>
501  static void mapToPhysicalFrame(ArrayPhysPoint & physPoints,
502  const ArrayRefPoint & refPoints,
503  const ArrayCell & cellWorkset,
504  const Teuchos::RCP< Basis< Scalar, FieldContainer<Scalar> > > HGRAD_Basis,
505  const int & whichCell = -1);
506 
562  /* template<class ArrayPhysPoint, class ArrayRefPoint, class ArrayCell>
563  static void mapToPhysicalFrameTemp(ArrayPhysPoint & physPoints,
564  const ArrayRefPoint & refPoints,
565  const ArrayCell & cellWorkset,
566  const shards::CellTopology & cellTopo,
567  const int & whichCell = -1);
568  */
569  template<class ArrayPhysPoint, class ArrayRefPoint, class ArrayCell, int refRank,int phyptsrank>
571 
630  template<class ArrayRefPoint, class ArrayPhysPoint, class ArrayCell>
631  static void mapToReferenceFrame(ArrayRefPoint & refPoints,
632  const ArrayPhysPoint & physPoints,
633  const ArrayCell & cellWorkset,
634  const shards::CellTopology & cellTopo,
635  const int & whichCell = -1);
636 
637 
638  template<class ArrayRefPoint, class ArrayPhysPoint, class ArrayCell>
639  static void mapToReferenceFrame(ArrayRefPoint & refPoints,
640  const ArrayPhysPoint & physPoints,
641  const ArrayCell & cellWorkset,
642  const Teuchos::RCP< Basis< Scalar, FieldContainer<Scalar> > > HGRAD_Basis,
643  const int & whichCell = -1);
644 
645 
691  template<class ArrayRefPoint, class ArrayInitGuess, class ArrayPhysPoint, class ArrayCell>
692  static void mapToReferenceFrameInitGuess(ArrayRefPoint & refPoints,
693  const ArrayInitGuess & initGuess,
694  const ArrayPhysPoint & physPoints,
695  const ArrayCell & cellWorkset,
696  const shards::CellTopology & cellTopo,
697  const int & whichCell = -1);
698 
699  template<class ArrayRefPoint, class ArrayInitGuess, class ArrayPhysPoint, class ArrayCell>
700  static void mapToReferenceFrameInitGuess(ArrayRefPoint & refPoints,
701  const ArrayInitGuess & initGuess,
702  const ArrayPhysPoint & physPoints,
703  const ArrayCell & cellWorkset,
704  const Teuchos::RCP< Basis< Scalar, FieldContainer<Scalar> > > HGRAD_Basis,
705  const int & whichCell = -1);
706 
707 
708 
759  template<class ArraySubcellPoint, class ArrayParamPoint>
760  static void mapToReferenceSubcell(ArraySubcellPoint & refSubcellPoints,
761  const ArrayParamPoint & paramPoints,
762  const int subcellDim,
763  const int subcellOrd,
764  const shards::CellTopology & parentCell);
765 
766 
767 
793  template<class ArrayEdgeTangent>
794  static void getReferenceEdgeTangent(ArrayEdgeTangent & refEdgeTangent,
795  const int & edgeOrd,
796  const shards::CellTopology & parentCell);
797 
798 
799 
836  template<class ArrayFaceTangentU, class ArrayFaceTangentV>
837  static void getReferenceFaceTangents(ArrayFaceTangentU & refFaceTanU,
838  ArrayFaceTangentV & refFaceTanV,
839  const int & faceOrd,
840  const shards::CellTopology & parentCell);
841 
842 
843 
906  template<class ArraySideNormal>
907  static void getReferenceSideNormal(ArraySideNormal & refSideNormal,
908  const int & sideOrd,
909  const shards::CellTopology & parentCell);
910 
911 
912 
951  template<class ArrayFaceNormal>
952  static void getReferenceFaceNormal(ArrayFaceNormal & refFaceNormal,
953  const int & faceOrd,
954  const shards::CellTopology & parentCell);
955 
956 
957 
987  template<class ArrayEdgeTangent, class ArrayJac>
988  static void getPhysicalEdgeTangents(ArrayEdgeTangent & edgeTangents,
989  const ArrayJac & worksetJacobians,
990  const int & worksetEdgeOrd,
991  const shards::CellTopology & parentCell);
992 
993  /* template<class ArrayEdgeTangent, class ArrayJac>
994  static void getPhysicalEdgeTangentsTemp(ArrayEdgeTangent & edgeTangents,
995  const ArrayJac & worksetJacobians,
996  const int & worksetEdgeOrd,
997  const shards::CellTopology & parentCell); */
998 
1038  template<class ArrayFaceTangentU, class ArrayFaceTangentV, class ArrayJac>
1039  static void getPhysicalFaceTangents(ArrayFaceTangentU & faceTanU,
1040  ArrayFaceTangentV & faceTanV,
1041  const ArrayJac & worksetJacobians,
1042  const int & worksetFaceOrd,
1043  const shards::CellTopology & parentCell);
1044 
1045  /* template<class ArrayFaceTangentU, class ArrayFaceTangentV, class ArrayJac>
1046  static void getPhysicalFaceTangentsTemp(ArrayFaceTangentU & faceTanU,
1047  ArrayFaceTangentV & faceTanV,
1048  const ArrayJac & worksetJacobians,
1049  const int & worksetFaceOrd,
1050  const shards::CellTopology & parentCell);
1051  */
1112  template<class ArraySideNormal, class ArrayJac>
1113  static void getPhysicalSideNormals(ArraySideNormal & sideNormals,
1114  const ArrayJac & worksetJacobians,
1115  const int & worksetSideOrd,
1116  const shards::CellTopology & parentCell);
1117 
1118 
1119 
1158  template<class ArrayFaceNormal, class ArrayJac>
1159  static void getPhysicalFaceNormals(ArrayFaceNormal & faceNormals,
1160  const ArrayJac & worksetJacobians,
1161  const int & worksetFaceOrd,
1162  const shards::CellTopology & parentCell);
1163 
1164  /* template<class ArrayFaceNormal, class ArrayJac>
1165  static void getPhysicalFaceNormalsTemp(ArrayFaceNormal & faceNormals,
1166  const ArrayJac & worksetJacobians,
1167  const int & worksetFaceOrd,
1168  const shards::CellTopology & parentCell);
1169  */
1170  //============================================================================================//
1171  // //
1172  // Inclusion tests //
1173  // //
1174  //============================================================================================//
1175 
1186  static int checkPointInclusion(const Scalar* point,
1187  const int pointDim,
1188  const shards::CellTopology & cellTopo,
1189  const double & threshold = INTREPID_THRESHOLD);
1190 
1191 
1192 
1205  template<class ArrayPoint>
1206  static int checkPointsetInclusion(const ArrayPoint & points,
1207  const shards::CellTopology & cellTopo,
1208  const double & threshold = INTREPID_THRESHOLD);
1209 
1210 
1211 
1238  template<class ArrayIncl, class ArrayPoint>
1239  static void checkPointwiseInclusion(ArrayIncl & inRefCell,
1240  const ArrayPoint & points,
1241  const shards::CellTopology & cellTopo,
1242  const double & threshold = INTREPID_THRESHOLD);
1243 
1244 
1245 
1281  template<class ArrayIncl, class ArrayPoint, class ArrayCell>
1282  static void checkPointwiseInclusion(ArrayIncl & inCell,
1283  const ArrayPoint & points,
1284  const ArrayCell & cellWorkset,
1285  const shards::CellTopology & cell,
1286  const int & whichCell = -1,
1287  const double & threshold = INTREPID_THRESHOLD);
1288 
1289 
1290 
1301  static const double* getReferenceVertex(const shards::CellTopology& cell,
1302  const int vertexOrd);
1303 
1304 
1305 
1320  template<class ArraySubcellVert>
1321  static void getReferenceSubcellVertices(ArraySubcellVert & subcellVertices,
1322  const int subcellDim,
1323  const int subcellOrd,
1324  const shards::CellTopology& parentCell);
1325 
1326 
1327 
1343  static const double* getReferenceNode(const shards::CellTopology& cell,
1344  const int nodeOrd);
1345 
1346 
1347 
1361  template<class ArraySubcellNode>
1362  static void getReferenceSubcellNodes(ArraySubcellNode& subcellNodes,
1363  const int subcellDim,
1364  const int subcellOrd,
1365  const shards::CellTopology& parentCell);
1366 
1367 
1368 
1374  static int hasReferenceCell(const shards::CellTopology & cellTopo);
1375 
1376 
1377 
1378  //============================================================================================//
1379  // //
1380  // Debug //
1381  // //
1382  //============================================================================================//
1383 
1384 
1390  static void printSubcellVertices(const int subcellDim,
1391  const int subcellOrd,
1392  const shards::CellTopology & parentCell);
1393 
1394 
1395 
1399  template<class ArrayCell>
1400  static void printWorksetSubcell(const ArrayCell & cellWorkset,
1401  const shards::CellTopology & parentCell,
1402  const int& pCellOrd,
1403  const int& subcellDim,
1404  const int& subcellOrd,
1405  const int& fieldWidth = 3);
1406 
1407  //============================================================================================//
1408  // //
1409  // Control Volume Coordinates //
1410  // //
1411  //============================================================================================//
1412 
1486  template<class ArrayCVCoord, class ArrayCellCoord>
1487  static void getSubCVCoords(ArrayCVCoord & subCVcoords, const ArrayCellCoord & cellCoords,
1488  const shards::CellTopology& primaryCell);
1489 
1496  template<class ArrayCent, class ArrayCellCoord>
1497  static void getBarycenter(ArrayCent & barycenter, const ArrayCellCoord & cellCoords);
1498 
1499 
1500  }; // class CellTools
1501 
1502 } // namespace Intrepid
1503 
1504 // include templated function definitions
1505 #ifdef INTREPID_OLD_KOKKOS_CODE
1506 #include <Intrepid_CellTools_Kokkos.hpp>
1507 #endif
1508 #include "Intrepid_CellToolsDef.hpp"
1509 
1510 #endif
1511 
1512 /***************************************************************************************************
1513  ** **
1514  ** D O C U M E N T A T I O N P A G E S **
1515  ** **
1516  **************************************************************************************************/
1517 
1798 #endif
static const double * getReferenceNode(const shards::CellTopology &cell, const int nodeOrd)
Retrieves the Cartesian coordinates of a reference cell node.
Header file for the Intrepid::HGRAD_WEDGE_I2_FEM class.
static void validateArguments_mapToReferenceFrame(const ArrayRefPoint &refPoints, const ArrayPhysPoint &physPoints, const ArrayCell &cellWorkset, const shards::CellTopology &cellTopo, const int &whichCell)
Validates arguments to Intrepid::CellTools::mapToReferenceFrame with default initial guess...
static void getReferenceEdgeTangent(ArrayEdgeTangent &refEdgeTangent, const int &edgeOrd, const shards::CellTopology &parentCell)
Computes constant tangent vectors to edges of 2D or 3D reference cells.
static void mapToReferenceSubcell(ArraySubcellPoint &refSubcellPoints, const ArrayParamPoint &paramPoints, const int subcellDim, const int subcellOrd, const shards::CellTopology &parentCell)
Computes parameterization maps of 1- and 2-subcells of reference cells.
Definition file for the Intrepid::CellTools class.
static void getPhysicalFaceNormals(ArrayFaceNormal &faceNormals, const ArrayJac &worksetJacobians, const int &worksetFaceOrd, const shards::CellTopology &parentCell)
Computes non-normalized normal vectors to physical faces in a face workset ; (see Subcell worksets fo...
static void setJacobianInv(ArrayJacInv &jacobianInv, const ArrayJac &jacobian)
Computes the inverse of the Jacobian matrix DF of the reference-to-physical frame map F...
static void getPhysicalEdgeTangents(ArrayEdgeTangent &edgeTangents, const ArrayJac &worksetJacobians, const int &worksetEdgeOrd, const shards::CellTopology &parentCell)
Computes non-normalized tangent vectors to physical edges in an edge workset ; (see Subcell worksets ...
static void mapToReferenceFrame(ArrayRefPoint &refPoints, const ArrayPhysPoint &physPoints, const ArrayCell &cellWorkset, const shards::CellTopology &cellTopo, const int &whichCell=-1)
Computes , the inverse of the reference-to-physical frame map using a default initial guess...
static void getReferenceSideNormal(ArraySideNormal &refSideNormal, const int &sideOrd, const shards::CellTopology &parentCell)
Computes constant normal vectors to sides of 2D or 3D reference cells.
static void getReferenceFaceTangents(ArrayFaceTangentU &refFaceTanU, ArrayFaceTangentV &refFaceTanV, const int &faceOrd, const shards::CellTopology &parentCell)
Computes pairs of constant tangent vectors to faces of a 3D reference cells.
static void printSubcellVertices(const int subcellDim, const int subcellOrd, const shards::CellTopology &parentCell)
Prints the reference space coordinates of the vertices of the specified subcell.
Header file for the Intrepid::HGRAD_TRI_C1_FEM class.
static void getReferenceFaceNormal(ArrayFaceNormal &refFaceNormal, const int &faceOrd, const shards::CellTopology &parentCell)
Computes constant normal vectors to faces of 3D reference cell.
Header file for utility class to provide multidimensional containers.
Contains definitions of custom data types in Intrepid.
static void mapToPhysicalFrame(ArrayPhysPoint &physPoints, const ArrayRefPoint &refPoints, const ArrayCell &cellWorkset, const shards::CellTopology &cellTopo, const int &whichCell=-1)
Computes F, the reference-to-physical frame map.
Header file for the Intrepid::HGRAD_TET_C2_FEM class.
static void getReferenceSubcellVertices(ArraySubcellVert &subcellVertices, const int subcellDim, const int subcellOrd, const shards::CellTopology &parentCell)
Retrieves the Cartesian coordinates of all vertices of a reference subcell.
static void getPhysicalSideNormals(ArraySideNormal &sideNormals, const ArrayJac &worksetJacobians, const int &worksetSideOrd, const shards::CellTopology &parentCell)
Computes non-normalized normal vectors to physical sides in a side workset .
static int checkPointsetInclusion(const ArrayPoint &points, const shards::CellTopology &cellTopo, const double &threshold=INTREPID_THRESHOLD)
Checks if a set of points belongs to a reference cell.
static int hasReferenceCell(const shards::CellTopology &cellTopo)
Checks if a cell topology has reference cell.
Header file for the Intrepid::HGRAD_TRI_C2_FEM class.
Intrepid utilities.
Header file for the abstract base class Intrepid::Basis.
Header file for the Intrepid::G_WEDGE_C2_FEM class.
static void mapToReferenceFrameInitGuess(ArrayRefPoint &refPoints, const ArrayInitGuess &initGuess, const ArrayPhysPoint &physPoints, const ArrayCell &cellWorkset, const shards::CellTopology &cellTopo, const int &whichCell=-1)
Computation of , the inverse of the reference-to-physical frame map using user-supplied initial guess...
static void setSubcellParametrization(FieldContainer< double > &subcellParametrization, const int subcellDim, const shards::CellTopology &parentCell)
Defines orientation-preserving parametrizations of reference edges and faces of cell topologies with ...
Header file for the Intrepid::HGRAD_HEX_I2_FEM class.
static void getReferenceSubcellNodes(ArraySubcellNode &subcellNodes, const int subcellDim, const int subcellOrd, const shards::CellTopology &parentCell)
Retrieves the Cartesian coordinates of all nodes of a reference subcell.
An abstract base class that defines interface for concrete basis implementations for Finite Element (...
static int checkPointInclusion(const Scalar *point, const int pointDim, const shards::CellTopology &cellTopo, const double &threshold=INTREPID_THRESHOLD)
Checks if a point belongs to a reference cell.
static const double * getReferenceVertex(const shards::CellTopology &cell, const int vertexOrd)
Retrieves the Cartesian coordinates of a reference cell vertex.
static void printWorksetSubcell(const ArrayCell &cellWorkset, const shards::CellTopology &parentCell, const int &pCellOrd, const int &subcellDim, const int &subcellOrd, const int &fieldWidth=3)
Prints the nodes of a subcell from a cell workset.
static void validateArguments_checkPointwiseInclusion(ArrayIncl &inCell, const ArrayPoint &physPoints, const ArrayCell &cellWorkset, const int &whichCell, const shards::CellTopology &cell)
Validates arguments to Intrepid::CellTools::checkPointwiseInclusion.
Computes the Jacobian matrix DF of the reference-to-physical frame map F.
static void validateArguments_setJacobianDetArgs(const ArrayJacDet &jacobianDet, const ArrayJac &jacobian)
Validates arguments to Intrepid::CellTools::setJacobianDet.
Header file for classes providing basic linear algebra functionality in 1D, 2D and 3D...
Header file for the Intrepid::HGRAD_HEX_C2_FEM class.
static void validateArguments_setJacobianInv(const ArrayJacInv &jacobianInv, const ArrayJac &jacobian)
Validates arguments to Intrepid::CellTools::setJacobianInv.
Computes F, the reference-to-physical frame map.
Header file for the Intrepid::HGRAD_TET_COMP12_FEM class.
static void validateArguments_mapToPhysicalFrame(const ArrayPhysPoint &physPoints, const ArrayRefPoint &refPoints, const ArrayCell &cellWorkset, const shards::CellTopology &cellTopo, const int &whichCell)
Validates arguments to Intrepid::CellTools::mapToPhysicalFrame.
static void checkPointwiseInclusion(ArrayIncl &inRefCell, const ArrayPoint &points, const shards::CellTopology &cellTopo, const double &threshold=INTREPID_THRESHOLD)
Checks every point in a set for inclusion in a reference cell.
static void validateArguments_setJacobian(const ArrayJac &jacobian, const ArrayPoint &points, const ArrayCell &cellWorkset, const int &whichCell, const shards::CellTopology &cellTopo)
Validates arguments to Intrepid::CellTools::setJacobian.
Header file for the Intrepid::HGRAD_PYR_I2_FEM class.
static const FieldContainer< double > & getSubcellParametrization(const int subcellDim, const shards::CellTopology &parentCell)
Returns array with the coefficients of the parametrization maps for the edges or faces of a reference...
static void getBarycenter(ArrayCent &barycenter, const ArrayCellCoord &cellCoords)
Compute cell barycenters.
static void getSubCVCoords(ArrayCVCoord &subCVcoords, const ArrayCellCoord &cellCoords, const shards::CellTopology &primaryCell)
Computes coordinates of sub-control volumes in each primary cell.
A stateless class for operations on cell data. Provides methods for:
static void setJacobianDet(ArrayJacDet &jacobianDet, const ArrayJac &jacobian)
Computes the determinant of the Jacobian matrix DF of the reference-to-physical frame map F...
CellTools()
Default constructor.
Header file for the Intrepid::HGRAD_LINE_C1_FEM class.
static void getPhysicalFaceTangents(ArrayFaceTangentU &faceTanU, ArrayFaceTangentV &faceTanV, const ArrayJac &worksetJacobians, const int &worksetFaceOrd, const shards::CellTopology &parentCell)
Computes non-normalized tangent vector pairs to physical faces in a face workset ; (see Subcell works...