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
48
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
70
75
78
82
84//#include "Intrepid2_HGRAD_WEDGE_I2_FEM.hpp"
85
86#include "Intrepid2_Data.hpp"
88
89namespace Intrepid2 {
90
91 //============================================================================================//
92 // //
93 // CellTools //
94 // //
95 //============================================================================================//
96
105
106
107 template<typename DeviceType>
108 class CellTools {
109 using ExecSpaceType = typename DeviceType::execution_space;
110 using MemSpaceType = typename DeviceType::memory_space;
111 public:
112
117 inline
118 static bool
119 hasReferenceCell( const shards::CellTopology cellTopo ) {
121 }
122
123 private:
124
128 template<typename outputValueType,
129 typename pointValueType>
130 static Teuchos::RCP<Basis<DeviceType,outputValueType,pointValueType> >
131 createHGradBasis( const shards::CellTopology cellTopo ) {
132 Teuchos::RCP<Basis<DeviceType,outputValueType,pointValueType> > r_val;
133
134 switch (cellTopo.getKey()) {
135 case shards::Line<2>::key: r_val = Teuchos::rcp(new Basis_HGRAD_LINE_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
136 case shards::Line<3>::key: r_val = Teuchos::rcp(new Basis_HGRAD_LINE_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
137 case shards::Triangle<3>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TRI_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
138 case shards::Quadrilateral<4>::key: r_val = Teuchos::rcp(new Basis_HGRAD_QUAD_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
139 case shards::Tetrahedron<4>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TET_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
140 case shards::Hexahedron<8>::key: r_val = Teuchos::rcp(new Basis_HGRAD_HEX_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
141 case shards::Wedge<6>::key: r_val = Teuchos::rcp(new Basis_HGRAD_WEDGE_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
142 case shards::Pyramid<5>::key: r_val = Teuchos::rcp(new Basis_HGRAD_PYR_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
143
144 case shards::Triangle<6>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TRI_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
145 case shards::Quadrilateral<9>::key: r_val = Teuchos::rcp(new Basis_HGRAD_QUAD_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
146 case shards::Tetrahedron<10>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TET_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
147 case shards::Tetrahedron<11>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TET_COMP12_FEM<DeviceType,outputValueType,pointValueType>()); break;
148 //case shards::Hexahedron<20>::key: r_val = Teuchos::rcp(new Basis_HGRAD_HEX_I2_FEM <DeviceType,outputValueType,pointValueType>()); break;
149 case shards::Hexahedron<27>::key: r_val = Teuchos::rcp(new Basis_HGRAD_HEX_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
150 //case shards::Wedge<15>::key: r_val = Teuchos::rcp(new Basis_HGRAD_WEDGE_I2_FEM <DeviceType,outputValueType,pointValueType>()); break;
151 case shards::Wedge<18>::key: r_val = Teuchos::rcp(new Basis_HGRAD_WEDGE_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
152 //case shards::Pyramid<13>::key: r_val = Teuchos::rcp(new Basis_HGRAD_PYR_I2_FEM <DeviceType,outputValueType,pointValueType>()); break;
153
154 case shards::Quadrilateral<8>::key:
155 case shards::Beam<2>::key:
156 case shards::Beam<3>::key:
157 case shards::ShellLine<2>::key:
158 case shards::ShellLine<3>::key:
159 case shards::ShellTriangle<3>::key:
160 case shards::ShellTriangle<6>::key:
161 case shards::ShellQuadrilateral<4>::key:
162 case shards::ShellQuadrilateral<8>::key:
163 case shards::ShellQuadrilateral<9>::key:
164 default: {
165 INTREPID2_TEST_FOR_EXCEPTION( true, std::invalid_argument,
166 ">>> ERROR (Intrepid2::CellTools::createHGradBasis): Cell topology not supported.");
167 }
168 }
169 return r_val;
170 }
171
172public:
173
176 CellTools() = default;
177
180 ~CellTools() = default;
181
182 //============================================================================================//
183 // //
184 // Jacobian, inverse Jacobian and Jacobian determinant //
185 // //
186 //============================================================================================//
187
223 template<typename jacobianValueType, class ...jacobianProperties,
224 typename pointValueType, class ...pointProperties,
225 typename WorksetType,
226 typename HGradBasisType>
227 static void
228 setJacobian( Kokkos::DynRankView<jacobianValueType,jacobianProperties...> jacobian,
229 const Kokkos::DynRankView<pointValueType,pointProperties...> points,
230 const WorksetType worksetCell,
231 const Teuchos::RCP<HGradBasisType> basis,
232 const int startCell=0, const int endCell=-1);
233
268 template<typename jacobianValueType, class ...jacobianProperties,
269 typename BasisGradientsType,
270 typename WorksetType>
271 static void
272 setJacobian( Kokkos::DynRankView<jacobianValueType,jacobianProperties...> jacobian,
273 const WorksetType worksetCell,
274 const BasisGradientsType gradients,
275 const int startCell=0, const int endCell=-1);
276
310
311 template<typename jacobianValueType, class ...jacobianProperties,
312 typename pointValueType, class ...pointProperties,
313 typename worksetCellValueType, class ...worksetCellProperties>
314 static void
315 setJacobian( Kokkos::DynRankView<jacobianValueType,jacobianProperties...> jacobian,
316 const Kokkos::DynRankView<pointValueType,pointProperties...> points,
317 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
318 const shards::CellTopology cellTopo ) {
319 using nonConstPointValueType = std::remove_const_t<pointValueType>;
321 setJacobian(jacobian,
322 points,
323 worksetCell,
324 basis);
325 }
326
337 template<typename jacobianInvValueType, class ...jacobianInvProperties,
338 typename jacobianValueType, class ...jacobianProperties>
339 static void
340 setJacobianInv( Kokkos::DynRankView<jacobianInvValueType,jacobianInvProperties...> jacobianInv,
341 const Kokkos::DynRankView<jacobianValueType,jacobianProperties...> jacobian );
342
353 template<typename jacobianDetValueType, class ...jacobianDetProperties,
354 typename jacobianValueType, class ...jacobianProperties>
355 static void
356 setJacobianDet( Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
357 const Kokkos::DynRankView<jacobianValueType,jacobianProperties...> jacobian );
358
364 template<class PointScalar>
366
372 template<class PointScalar>
374
380 template<class PointScalar>
381 static void setJacobianDet( Data<PointScalar,DeviceType> & jacobianDet,
382 const Data<PointScalar,DeviceType> & jacobian);
383
389 template<class PointScalar>
390 static void setJacobianDetInv( Data<PointScalar,DeviceType> & jacobianDet,
391 const Data<PointScalar,DeviceType> & jacobian);
392
398 template<class PointScalar>
399 static void setJacobianInv( Data<PointScalar,DeviceType> & jacobianInv,
400 const Data<PointScalar,DeviceType> & jacobian);
401
408 template<class PointScalar>
409 static void setJacobianDividedByDet( Data<PointScalar,DeviceType> & jacobianDividedByDet,
410 const Data<PointScalar,DeviceType> & jacobian,
411 const Data<PointScalar,DeviceType> & jacobianDetInv);
412
413 //============================================================================================//
414 // //
415 // Node information //
416 // //
417 //============================================================================================//
418
419 // the node information can be used inside of kokkos functor and needs kokkos inline and
420 // exception should be an abort. for now, let's not decorate
421
429 template<typename cellCenterValueType, class ...cellCenterProperties>
430 static void
431 getReferenceCellCenter( Kokkos::DynRankView<cellCenterValueType,cellCenterProperties...> cellCenter,
432 const shards::CellTopology cell );
433
442 template<typename cellVertexValueType, class ...cellVertexProperties>
443 static void
444 getReferenceVertex( Kokkos::DynRankView<cellVertexValueType,cellVertexProperties...> cellVertex,
445 const shards::CellTopology cell,
446 const ordinal_type vertexOrd );
447
448
463 template<typename subcellVertexValueType, class ...subcellVertexProperties>
464 static void
465 getReferenceSubcellVertices( Kokkos::DynRankView<subcellVertexValueType,subcellVertexProperties...> subcellVertices,
466 const ordinal_type subcellDim,
467 const ordinal_type subcellOrd,
468 const shards::CellTopology parentCell );
469
470
471
487 template<typename cellNodeValueType, class ...cellNodeProperties>
488 static void
489 getReferenceNode( Kokkos::DynRankView<cellNodeValueType,cellNodeProperties...> cellNode,
490 const shards::CellTopology cell,
491 const ordinal_type nodeOrd );
492
493
494
508 template<typename subcellNodeValueType, class ...subcellNodeProperties>
509 static void
510 getReferenceSubcellNodes( Kokkos::DynRankView<subcellNodeValueType,subcellNodeProperties...> subcellNodes,
511 const ordinal_type subcellDim,
512 const ordinal_type subcellOrd,
513 const shards::CellTopology parentCell );
514
540 template<typename refEdgeTangentValueType, class ...refEdgeTangentProperties>
541 static void
542 getReferenceEdgeTangent( Kokkos::DynRankView<refEdgeTangentValueType,refEdgeTangentProperties...> refEdgeTangent,
543 const ordinal_type edgeOrd,
544 const shards::CellTopology parentCell );
545
582 template<typename refFaceTanValueType, class ...refFaceTanProperties>
583 static void
584 getReferenceFaceTangents( Kokkos::DynRankView<refFaceTanValueType,refFaceTanProperties...> refFaceTanU,
585 Kokkos::DynRankView<refFaceTanValueType,refFaceTanProperties...> refFaceTanV,
586 const ordinal_type faceOrd,
587 const shards::CellTopology parentCell );
588
651 template<typename refSideNormalValueType, class ...refSideNormalProperties>
652 static void
653 getReferenceSideNormal( Kokkos::DynRankView<refSideNormalValueType,refSideNormalProperties...> refSideNormal,
654 const ordinal_type sideOrd,
655 const shards::CellTopology parentCell );
656
695 template<typename refFaceNormalValueType, class ...refFaceNormalProperties>
696 static void
697 getReferenceFaceNormal( Kokkos::DynRankView<refFaceNormalValueType,refFaceNormalProperties...> refFaceNormal,
698 const ordinal_type faceOrd,
699 const shards::CellTopology parentCell );
700
730 template<typename edgeTangentValueType, class ...edgeTangentProperties,
731 typename worksetJacobianValueType, class ...worksetJacobianProperties>
732 static void
733 getPhysicalEdgeTangents( Kokkos::DynRankView<edgeTangentValueType,edgeTangentProperties...> edgeTangents,
734 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
735 const ordinal_type worksetEdgeOrd,
736 const shards::CellTopology parentCell );
737
738
750 template<typename edgeTangentValueType, class ...edgeTangentProperties,
751 typename worksetJacobianValueType, class ...worksetJacobianProperties,
752 typename edgeOrdValueType, class ...edgeOrdProperties>
753 static void
754 getPhysicalEdgeTangents( Kokkos::DynRankView<edgeTangentValueType,edgeTangentProperties...> edgeTangents,
755 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
756 const Kokkos::DynRankView<edgeOrdValueType,edgeOrdProperties...> worksetEdgeOrds,
757 const shards::CellTopology parentCell );
758
798 template<typename faceTanValueType, class ...faceTanProperties,
799 typename worksetJacobianValueType, class ...worksetJacobianProperties>
800 static void
801 getPhysicalFaceTangents( Kokkos::DynRankView<faceTanValueType,faceTanProperties...> faceTanU,
802 Kokkos::DynRankView<faceTanValueType,faceTanProperties...> faceTanV,
803 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
804 const ordinal_type worksetFaceOrd,
805 const shards::CellTopology parentCell );
806
807
820 template<typename faceTanValueType, class ...faceTanProperties,
821 typename worksetJacobianValueType, class ...worksetJacobianProperties,
822 typename faceOrdValueType, class ...faceOrdProperties>
823 static void
824 getPhysicalFaceTangents( Kokkos::DynRankView<faceTanValueType,faceTanProperties...> faceTanU,
825 Kokkos::DynRankView<faceTanValueType,faceTanProperties...> faceTanV,
826 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
827 const Kokkos::DynRankView<faceOrdValueType,faceOrdProperties...> worksetFaceOrds,
828 const shards::CellTopology parentCell );
829
830
831
892 template<typename sideNormalValueType, class ...sideNormalProperties,
893 typename worksetJacobianValueType, class ...worksetJacobianProperties>
894 static void
895 getPhysicalSideNormals( Kokkos::DynRankView<sideNormalValueType,sideNormalProperties...> sideNormals,
896 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
897 const ordinal_type worksetSideOrd,
898 const shards::CellTopology parentCell );
899
900
911 template<typename sideNormalValueType, class ...sideNormalProperties,
912 typename worksetJacobianValueType, class ...worksetJacobianProperties,
913 typename edgeOrdValueType, class ...edgeOrdProperties>
914 static void
915 getPhysicalSideNormals( Kokkos::DynRankView<sideNormalValueType,sideNormalProperties...> sideNormals,
916 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
917 const Kokkos::DynRankView<edgeOrdValueType,edgeOrdProperties...> worksetSideOrds,
918 const shards::CellTopology parentCell );
919
958 template<typename faceNormalValueType, class ...faceNormalProperties,
959 typename worksetJacobianValueType, class ...worksetJacobianProperties>
960 static void
961 getPhysicalFaceNormals( Kokkos::DynRankView<faceNormalValueType,faceNormalProperties...> faceNormals,
962 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
963 const ordinal_type worksetFaceOrd,
964 const shards::CellTopology parentCell );
965
966
977
978 template<typename faceNormalValueType, class ...faceNormalProperties,
979 typename worksetJacobianValueType, class ...worksetJacobianProperties,
980 typename faceOrdValueType, class ...faceOrdProperties>
981 static void
982 getPhysicalFaceNormals( Kokkos::DynRankView<faceNormalValueType,faceNormalProperties...> faceNormals,
983 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
984 const Kokkos::DynRankView<faceOrdValueType,faceOrdProperties...> worksetFaceOrds,
985 const shards::CellTopology parentCell );
986
987 //============================================================================================//
988 // //
989 // Reference-to-physical frame mapping and its inverse //
990 // //
991 //============================================================================================//
992
1029 template<typename physPointValueType, class ...physPointProperties,
1030 typename refPointValueType, class ...refPointProperties,
1031 typename WorksetType,
1032 typename HGradBasisPtrType>
1033 static void
1034 mapToPhysicalFrame( Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1035 const Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
1036 const WorksetType worksetCell,
1037 const HGradBasisPtrType basis );
1038
1079 template<typename physPointValueType, class ...physPointProperties,
1080 typename refPointValueType, class ...refPointProperties,
1081 typename worksetCellValueType, class ...worksetCellProperties>
1082 static void
1083 mapToPhysicalFrame( Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1084 const Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
1085 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1086 const shards::CellTopology cellTopo ) {
1087 using nonConstRefPointValueType = std::remove_const_t<refPointValueType>;
1089 mapToPhysicalFrame(physPoints,
1090 refPoints,
1091 worksetCell,
1092 basis);
1093 }
1094
1095
1146 template<typename refSubcellPointValueType, class ...refSubcellPointProperties,
1147 typename paramPointValueType, class ...paramPointProperties>
1148 static void
1149 mapToReferenceSubcell( Kokkos::DynRankView<refSubcellPointValueType,refSubcellPointProperties...> refSubcellPoints,
1150 const Kokkos::DynRankView<paramPointValueType,paramPointProperties...> paramPoints,
1151 const ordinal_type subcellDim,
1152 const ordinal_type subcellOrd,
1153 const shards::CellTopology parentCell );
1154
1155
1156 //============================================================================================//
1157 // //
1158 // Physical-to-reference frame mapping and its inverse //
1159 // //
1160 //============================================================================================//
1161
1162
1206 template<typename refPointValueType, class ...refPointProperties,
1207 typename physPointValueType, class ...physPointProperties,
1208 typename worksetCellValueType, class ...worksetCellProperties>
1209 static void
1210 mapToReferenceFrame( Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
1211 const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1212 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1213 const shards::CellTopology cellTopo );
1214
1241 template<typename refPointValueType, class ...refPointProperties,
1242 typename initGuessValueType, class ...initGuessProperties,
1243 typename physPointValueType, class ...physPointProperties,
1244 typename worksetCellValueType, class ...worksetCellProperties,
1245 typename HGradBasisPtrType>
1246 static void
1247 mapToReferenceFrameInitGuess( Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
1248 const Kokkos::DynRankView<initGuessValueType,initGuessProperties...> initGuess,
1249 const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1250 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1251 const HGradBasisPtrType basis );
1252
1253
1284 template<typename refPointValueType, class ...refPointProperties,
1285 typename initGuessValueType, class ...initGuessProperties,
1286 typename physPointValueType, class ...physPointProperties,
1287 typename worksetCellValueType, class ...worksetCellProperties>
1288 static void
1289 mapToReferenceFrameInitGuess( Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
1290 const Kokkos::DynRankView<initGuessValueType,initGuessProperties...> initGuess,
1291 const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1292 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1293 const shards::CellTopology cellTopo ) {
1296 initGuess,
1297 physPoints,
1298 worksetCell,
1299 basis);
1300 }
1301
1302
1303 //============================================================================================//
1304 // //
1305 // Control Volume Coordinates //
1306 // //
1307 //============================================================================================//
1308
1382 template<typename subcvCoordValueType, class ...subcvCoordProperties,
1383 typename cellCoordValueType, class ...cellCoordProperties>
1384 static void
1385 getSubcvCoords( Kokkos::DynRankView<subcvCoordValueType,subcvCoordProperties...> subcvCoords,
1386 const Kokkos::DynRankView<cellCoordValueType,cellCoordProperties...> cellCoords,
1387 const shards::CellTopology primaryCell );
1388
1389 //============================================================================================//
1390 // //
1391 // Inclusion tests //
1392 // //
1393 //============================================================================================//
1394
1404 template<typename pointValueType, class ...pointProperties>
1405 static bool
1406 checkPointInclusion( const Kokkos::DynRankView<pointValueType,pointProperties...> point,
1407 const shards::CellTopology cellTopo,
1408 const double thres = threshold() );
1409
1410 // template<class ArrayPoint>
1411 // static int checkPointsetInclusion(const ArrayPoint & points,
1412 // const shards::CellTopology & cellTopo,
1413 // const double & threshold = INTREPID2_THRESHOLD);
1414
1415
1416
1443 template<typename inCellValueType, class ...inCellProperties,
1444 typename pointValueType, class ...pointProperties>
1445 static void checkPointwiseInclusion( Kokkos::DynRankView<inCellValueType,inCellProperties...> inCell,
1446 const Kokkos::DynRankView<pointValueType,pointProperties...> points,
1447 const shards::CellTopology cellTopo,
1448 const double thres = threshold() );
1449
1471 template<typename inCellValueType, class ...inCellProperties,
1472 typename pointValueType, class ...pointProperties,
1473 typename cellWorksetValueType, class ...cellWorksetProperties>
1474 static void checkPointwiseInclusion( Kokkos::DynRankView<inCellValueType,inCellProperties...> inCell,
1475 const Kokkos::DynRankView<pointValueType,pointProperties...> points,
1476 const Kokkos::DynRankView<cellWorksetValueType,cellWorksetProperties...> cellWorkset,
1477 const shards::CellTopology cellTopo,
1478 const double thres = threshold() );
1479
1480
1481 // //============================================================================================//
1482 // // //
1483 // // Debug //
1484 // // //
1485 // //============================================================================================//
1486
1487
1488 // /** \brief Prints the reference space coordinates of the vertices of the specified subcell
1489 // \param subcellDim [in] - dimension of the subcell where points are mapped to
1490 // \param subcellOrd [in] - subcell ordinal
1491 // \param parentCell [in] - cell topology of the parent cell.
1492 // */
1493 // static void printSubcellVertices(const int subcellDim,
1494 // const int subcellOrd,
1495 // const shards::CellTopology & parentCell);
1496
1497
1498
1499 // /** \brief Prints the nodes of a subcell from a cell workset
1500
1501 // */
1502 // template<class ArrayCell>
1503 // static void printWorksetSubcell(const ArrayCell & cellWorkset,
1504 // const shards::CellTopology & parentCell,
1505 // const int& pCellOrd,
1506 // const int& subcellDim,
1507 // const int& subcellOrd,
1508 // const int& fieldWidth = 3);
1509 };
1510
1511 //============================================================================================//
1512 // //
1513 // Validation of input/output arguments for CellTools methods //
1514 // //
1515 //============================================================================================//
1516
1523 template<typename jacobianViewType,
1524 typename PointViewType,
1525 typename worksetCellViewType>
1526 static void
1527 CellTools_setJacobianArgs( const jacobianViewType jacobian,
1528 const PointViewType points,
1529 const worksetCellViewType worksetCell,
1530 const shards::CellTopology cellTopo );
1531
1536 template<typename jacobianInvViewType,
1537 typename jacobianViewType>
1538 static void
1539 CellTools_setJacobianInvArgs( const jacobianInvViewType jacobianInv,
1540 const jacobianViewType jacobian );
1541
1542
1547 template<typename jacobianDetViewType,
1548 typename jacobianViewType>
1549 static void
1550 CellTools_setJacobianDetArgs( const jacobianDetViewType jacobianDet,
1551 const jacobianViewType jacobian );
1552
1553
1560 template<typename physPointViewType,
1561 typename refPointViewType,
1562 typename worksetCellViewType>
1563 static void
1564 CellTools_mapToPhysicalFrameArgs( const physPointViewType physPoints,
1565 const refPointViewType refPoints,
1566 const worksetCellViewType worksetCell,
1567 const shards::CellTopology cellTopo );
1568
1569
1576 template<typename refPointViewType,
1577 typename physPointViewType,
1578 typename worksetCellViewType>
1579 static void
1580 CellTools_mapToReferenceFrameArgs( const refPointViewType refPoints,
1581 const physPointViewType physPoints,
1582 const worksetCellViewType worksetCell,
1583 const shards::CellTopology cellTopo );
1584
1585
1586
1594 template<typename refPointViewType,
1595 typename initGuessViewType,
1596 typename physPointViewType,
1597 typename worksetCellViewType>
1598 static void
1599 CellTools_mapToReferenceFrameInitGuess( const refPointViewType refPoints,
1600 const initGuessViewType initGuess,
1601 const physPointViewType physPoints,
1602 const worksetCellViewType worksetCell,
1603 const shards::CellTopology cellTopo );
1604
1605 // /** \brief Validates arguments to Intrepid2::CellTools::checkPointwiseInclusion
1606 // \param inCell [out] - rank-1 (P) array required
1607 // \param physPoints [in] - rank-2 (P,D) array required
1608 // \param cellWorkset [in] - rank-3 (C,N,D) array required
1609 // \param whichCell [in] - 0 <= whichCell < C required
1610 // \param cellTopo [in] - cell topology with a reference cell required
1611 // */
1612 // template<class ArrayIncl, class ArrayPoint, class ArrayCell>
1613 // static void
1614 // validateArguments_checkPointwiseInclusion(ArrayIncl & inCell,
1615 // const ArrayPoint & physPoints,
1616 // const ArrayCell & cellWorkset,
1617 // const int & whichCell,
1618 // const shards::CellTopology & cell);
1619
1620
1621
1622}
1623
1625
1628
1632
1634
1635// not yet converted ...
1637
1638// #include "Intrepid2_CellToolsDefDebug.hpp"
1639
1640
1641#endif
1642
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.
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.
static void CellTools_mapToPhysicalFrameArgs(const physPointViewType physPoints, const refPointViewType refPoints, const worksetCellViewType worksetCell, const shards::CellTopology cellTopo)
Validates arguments to Intrepid2::CellTools::mapToPhysicalFrame.
static void CellTools_setJacobianInvArgs(const jacobianInvViewType jacobianInv, const jacobianViewType jacobian)
Validates arguments to Intrepid2::CellTools::setJacobianInv.
static void CellTools_setJacobianDetArgs(const jacobianDetViewType jacobianDet, const jacobianViewType jacobian)
Validates arguments to Intrepid2::CellTools::setJacobianDet.
static 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.
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_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.
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 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 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 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(Kokkos::DynRankView< physPointValueType, physPointProperties... > physPoints, const Kokkos::DynRankView< refPointValueType, refPointProperties... > refPoints, const WorksetType worksetCell, const HGradBasisPtrType basis)
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 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 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 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 setJacobianDetInv(Data< PointScalar, DeviceType > &jacobianDet, const Data< PointScalar, DeviceType > &jacobian)
Computes reciprocals of determinants corresponding to the Jacobians in the Data container provided.
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 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 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.
static bool hasReferenceCell(const shards::CellTopology cellTopo)
Checks if a cell topology has reference cell.
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.
static void setJacobian(Kokkos::DynRankView< jacobianValueType, jacobianProperties... > jacobian, const Kokkos::DynRankView< pointValueType, pointProperties... > 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.
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...
CellTools()=default
Default constructor.
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 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.
~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 getReferenceCellCenter(Kokkos::DynRankView< cellCenterValueType, cellCenterProperties... > cellCenter, const shards::CellTopology cell)
Computes the Cartesian coordinates of reference cell barycenter.
static void getReferenceFaceTangents(Kokkos::DynRankView< refFaceTanValueType, refFaceTanProperties... > refFaceTanU, Kokkos::DynRankView< refFaceTanValueType, refFaceTanProperties... > refFaceTanV, const ordinal_type faceOrd, const shards::CellTopology parentCell)
Computes pairs of constant tangent vectors to faces of a 3D reference cells.
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 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 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 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.