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