Intrepid2
Intrepid2_CellTools.hpp
Go to the documentation of this file.
1// @HEADER
2// ************************************************************************
3//
4// Intrepid2 Package
5// Copyright (2007) Sandia Corporation
6//
7// Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
8// license for use of this work by or on behalf of the U.S. Government.
9//
10// Redistribution and use in source and binary forms, with or without
11// modification, are permitted provided that the following conditions are
12// met:
13//
14// 1. Redistributions of source code must retain the above copyright
15// notice, this list of conditions and the following disclaimer.
16//
17// 2. Redistributions in binary form must reproduce the above copyright
18// notice, this list of conditions and the following disclaimer in the
19// documentation and/or other materials provided with the distribution.
20//
21// 3. Neither the name of the Corporation nor the names of the
22// contributors may be used to endorse or promote products derived from
23// this software without specific prior written permission.
24//
25// THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
26// EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
27// IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
28// PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
29// CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
30// EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
31// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
32// PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
33// LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
34// NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
35// SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
36//
37// Questions? Contact Kyungjoo Kim (kyukim@sandia.gov), or
38// Mauro Perego (mperego@sandia.gov)
39//
40// ************************************************************************
41// @HEADER
42
49#ifndef __INTREPID2_CELLTOOLS_HPP__
50#define __INTREPID2_CELLTOOLS_HPP__
51
52#include "Intrepid2_ConfigDefs.hpp"
53
54#include "Shards_CellTopology.hpp"
55#include "Shards_BasicTopologies.hpp"
56
57#include "Teuchos_RCP.hpp"
58
59#include "Intrepid2_Types.hpp"
60#include "Intrepid2_Utils.hpp"
61
63
64#include "Intrepid2_Basis.hpp"
65
69
74
77
81
83//#include "Intrepid2_HGRAD_WEDGE_I2_FEM.hpp"
84
85#include "Intrepid2_Data.hpp"
87
88namespace Intrepid2 {
89
90 //============================================================================================//
91 // //
92 // CellTools //
93 // //
94 //============================================================================================//
95
106 template<typename DeviceType>
107 class CellTools {
108 using ExecSpaceType = typename DeviceType::execution_space;
109 using MemSpaceType = typename DeviceType::memory_space;
110 public:
111
116 inline
117 static bool
118 hasReferenceCell( const shards::CellTopology cellTopo ) {
120 }
121
122 private:
123
127 template<typename outputValueType,
128 typename pointValueType>
129 static Teuchos::RCP<Basis<DeviceType,outputValueType,pointValueType> >
130 createHGradBasis( const shards::CellTopology cellTopo ) {
131 Teuchos::RCP<Basis<DeviceType,outputValueType,pointValueType> > r_val;
132
133 switch (cellTopo.getKey()) {
134 case shards::Line<2>::key: r_val = Teuchos::rcp(new Basis_HGRAD_LINE_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
135 case shards::Triangle<3>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TRI_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
136 case shards::Quadrilateral<4>::key: r_val = Teuchos::rcp(new Basis_HGRAD_QUAD_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
137 case shards::Tetrahedron<4>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TET_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
138 case shards::Hexahedron<8>::key: r_val = Teuchos::rcp(new Basis_HGRAD_HEX_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
139 case shards::Wedge<6>::key: r_val = Teuchos::rcp(new Basis_HGRAD_WEDGE_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
140 case shards::Pyramid<5>::key: r_val = Teuchos::rcp(new Basis_HGRAD_PYR_C1_FEM <DeviceType,outputValueType,pointValueType>()); break;
141
142 case shards::Triangle<6>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TRI_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
143 case shards::Quadrilateral<9>::key: r_val = Teuchos::rcp(new Basis_HGRAD_QUAD_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
144 case shards::Tetrahedron<10>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TET_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
145 case shards::Tetrahedron<11>::key: r_val = Teuchos::rcp(new Basis_HGRAD_TET_COMP12_FEM<DeviceType,outputValueType,pointValueType>()); break;
146 //case shards::Hexahedron<20>::key: r_val = Teuchos::rcp(new Basis_HGRAD_HEX_I2_FEM <DeviceType,outputValueType,pointValueType>()); break;
147 case shards::Hexahedron<27>::key: r_val = Teuchos::rcp(new Basis_HGRAD_HEX_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
148 //case shards::Wedge<15>::key: r_val = Teuchos::rcp(new Basis_HGRAD_WEDGE_I2_FEM <DeviceType,outputValueType,pointValueType>()); break;
149 case shards::Wedge<18>::key: r_val = Teuchos::rcp(new Basis_HGRAD_WEDGE_C2_FEM <DeviceType,outputValueType,pointValueType>()); break;
150 //case shards::Pyramid<13>::key: r_val = Teuchos::rcp(new Basis_HGRAD_PYR_I2_FEM <DeviceType,outputValueType,pointValueType>()); break;
151
152 case shards::Quadrilateral<8>::key:
153 case shards::Line<3>::key:
154 case shards::Beam<2>::key:
155 case shards::Beam<3>::key:
156 case shards::ShellLine<2>::key:
157 case shards::ShellLine<3>::key:
158 case shards::ShellTriangle<3>::key:
159 case shards::ShellTriangle<6>::key:
160 case shards::ShellQuadrilateral<4>::key:
161 case shards::ShellQuadrilateral<8>::key:
162 case shards::ShellQuadrilateral<9>::key:
163 default: {
164 INTREPID2_TEST_FOR_EXCEPTION( true, std::invalid_argument,
165 ">>> ERROR (Intrepid2::CellTools::createHGradBasis): Cell topology not supported.");
166 }
167 }
168 return r_val;
169 }
170
171public:
172
175 CellTools() = default;
176
179 ~CellTools() = default;
180
181 //============================================================================================//
182 // //
183 // Jacobian, inverse Jacobian and Jacobian determinant //
184 // //
185 //============================================================================================//
186
222 template<typename jacobianValueType, class ...jacobianProperties,
223 typename pointValueType, class ...pointProperties,
224 typename WorksetType,
225 typename HGradBasisType>
226 static void
227 setJacobian( Kokkos::DynRankView<jacobianValueType,jacobianProperties...> jacobian,
228 const Kokkos::DynRankView<pointValueType,pointProperties...> points,
229 const WorksetType worksetCell,
230 const Teuchos::RCP<HGradBasisType> basis,
231 const int startCell=0, const int endCell=-1);
232
267 template<typename jacobianValueType, class ...jacobianProperties,
268 typename BasisGradientsType,
269 typename WorksetType>
270 static void
271 setJacobian( Kokkos::DynRankView<jacobianValueType,jacobianProperties...> jacobian,
272 const WorksetType worksetCell,
273 const BasisGradientsType gradients,
274 const int startCell=0, const int endCell=-1);
275
310 template<typename jacobianValueType, class ...jacobianProperties,
311 typename pointValueType, class ...pointProperties,
312 typename worksetCellValueType, class ...worksetCellProperties>
313 static void
314 setJacobian( Kokkos::DynRankView<jacobianValueType,jacobianProperties...> jacobian,
315 const Kokkos::DynRankView<pointValueType,pointProperties...> points,
316 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
317 const shards::CellTopology cellTopo ) {
318 auto basis = createHGradBasis<pointValueType,pointValueType>(cellTopo);
319 setJacobian(jacobian,
320 points,
321 worksetCell,
322 basis);
323 }
324
335 template<typename jacobianInvValueType, class ...jacobianInvProperties,
336 typename jacobianValueType, class ...jacobianProperties>
337 static void
338 setJacobianInv( Kokkos::DynRankView<jacobianInvValueType,jacobianInvProperties...> jacobianInv,
339 const Kokkos::DynRankView<jacobianValueType,jacobianProperties...> jacobian );
340
351 template<typename jacobianDetValueType, class ...jacobianDetProperties,
352 typename jacobianValueType, class ...jacobianProperties>
353 static void
354 setJacobianDet( Kokkos::DynRankView<jacobianDetValueType,jacobianDetProperties...> jacobianDet,
355 const Kokkos::DynRankView<jacobianValueType,jacobianProperties...> jacobian );
356
362 template<class PointScalar>
364
370 template<class PointScalar>
372
378 template<class PointScalar>
379 static void setJacobianDet( Data<PointScalar,DeviceType> & jacobianDet,
380 const Data<PointScalar,DeviceType> & jacobian);
381
387 template<class PointScalar>
388 static void setJacobianInv( Data<PointScalar,DeviceType> & jacobianInv,
389 const Data<PointScalar,DeviceType> & jacobian);
390
391 //============================================================================================//
392 // //
393 // Node information //
394 // //
395 //============================================================================================//
396
397 // the node information can be used inside of kokkos functor and needs kokkos inline and
398 // exception should be an abort. for now, let's not decorate
399
407 template<typename cellCenterValueType, class ...cellCenterProperties>
408 static void
409 getReferenceCellCenter( Kokkos::DynRankView<cellCenterValueType,cellCenterProperties...> cellCenter,
410 const shards::CellTopology cell );
411
420 template<typename cellVertexValueType, class ...cellVertexProperties>
421 static void
422 getReferenceVertex( Kokkos::DynRankView<cellVertexValueType,cellVertexProperties...> cellVertex,
423 const shards::CellTopology cell,
424 const ordinal_type vertexOrd );
425
426
441 template<typename subcellVertexValueType, class ...subcellVertexProperties>
442 static void
443 getReferenceSubcellVertices( Kokkos::DynRankView<subcellVertexValueType,subcellVertexProperties...> subcellVertices,
444 const ordinal_type subcellDim,
445 const ordinal_type subcellOrd,
446 const shards::CellTopology parentCell );
447
448
449
465 template<typename cellNodeValueType, class ...cellNodeProperties>
466 static void
467 getReferenceNode( Kokkos::DynRankView<cellNodeValueType,cellNodeProperties...> cellNode,
468 const shards::CellTopology cell,
469 const ordinal_type nodeOrd );
470
471
472
486 template<typename subcellNodeValueType, class ...subcellNodeProperties>
487 static void
488 getReferenceSubcellNodes( Kokkos::DynRankView<subcellNodeValueType,subcellNodeProperties...> subcellNodes,
489 const ordinal_type subcellDim,
490 const ordinal_type subcellOrd,
491 const shards::CellTopology parentCell );
492
518 template<typename refEdgeTangentValueType, class ...refEdgeTangentProperties>
519 static void
520 getReferenceEdgeTangent( Kokkos::DynRankView<refEdgeTangentValueType,refEdgeTangentProperties...> refEdgeTangent,
521 const ordinal_type edgeOrd,
522 const shards::CellTopology parentCell );
523
560 template<typename refFaceTanValueType, class ...refFaceTanProperties>
561 static void
562 getReferenceFaceTangents( Kokkos::DynRankView<refFaceTanValueType,refFaceTanProperties...> refFaceTanU,
563 Kokkos::DynRankView<refFaceTanValueType,refFaceTanProperties...> refFaceTanV,
564 const ordinal_type faceOrd,
565 const shards::CellTopology parentCell );
566
629 template<typename refSideNormalValueType, class ...refSideNormalProperties>
630 static void
631 getReferenceSideNormal( Kokkos::DynRankView<refSideNormalValueType,refSideNormalProperties...> refSideNormal,
632 const ordinal_type sideOrd,
633 const shards::CellTopology parentCell );
634
673 template<typename refFaceNormalValueType, class ...refFaceNormalProperties>
674 static void
675 getReferenceFaceNormal( Kokkos::DynRankView<refFaceNormalValueType,refFaceNormalProperties...> refFaceNormal,
676 const ordinal_type faceOrd,
677 const shards::CellTopology parentCell );
678
708 template<typename edgeTangentValueType, class ...edgeTangentProperties,
709 typename worksetJacobianValueType, class ...worksetJacobianProperties>
710 static void
711 getPhysicalEdgeTangents( Kokkos::DynRankView<edgeTangentValueType,edgeTangentProperties...> edgeTangents,
712 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
713 const ordinal_type worksetEdgeOrd,
714 const shards::CellTopology parentCell );
715
755 template<typename faceTanValueType, class ...faceTanProperties,
756 typename worksetJacobianValueType, class ...worksetJacobianProperties>
757 static void
758 getPhysicalFaceTangents( Kokkos::DynRankView<faceTanValueType,faceTanProperties...> faceTanU,
759 Kokkos::DynRankView<faceTanValueType,faceTanProperties...> faceTanV,
760 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
761 const ordinal_type worksetFaceOrd,
762 const shards::CellTopology parentCell );
763
824 template<typename sideNormalValueType, class ...sideNormalProperties,
825 typename worksetJacobianValueType, class ...worksetJacobianProperties>
826 static void
827 getPhysicalSideNormals( Kokkos::DynRankView<sideNormalValueType,sideNormalProperties...> sideNormals,
828 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
829 const ordinal_type worksetSideOrd,
830 const shards::CellTopology parentCell );
831
870 template<typename faceNormalValueType, class ...faceNormalProperties,
871 typename worksetJacobianValueType, class ...worksetJacobianProperties>
872 static void
873 getPhysicalFaceNormals( Kokkos::DynRankView<faceNormalValueType,faceNormalProperties...> faceNormals,
874 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
875 const ordinal_type worksetFaceOrd,
876 const shards::CellTopology parentCell );
877
878 //============================================================================================//
879 // //
880 // Reference-to-physical frame mapping and its inverse //
881 // //
882 //============================================================================================//
883
920 template<typename physPointValueType, class ...physPointProperties,
921 typename refPointValueType, class ...refPointProperties,
922 typename WorksetType,
923 typename HGradBasisPtrType>
924 static void
925 mapToPhysicalFrame( Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
926 const Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
927 const WorksetType worksetCell,
928 const HGradBasisPtrType basis );
929
970 template<typename physPointValueType, class ...physPointProperties,
971 typename refPointValueType, class ...refPointProperties,
972 typename worksetCellValueType, class ...worksetCellProperties>
973 static void
974 mapToPhysicalFrame( Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
975 const Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
976 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
977 const shards::CellTopology cellTopo ) {
978 auto basis = createHGradBasis<refPointValueType,refPointValueType>(cellTopo);
979 mapToPhysicalFrame(physPoints,
980 refPoints,
981 worksetCell,
982 basis);
983 }
984
985
1036 template<typename refSubcellPointValueType, class ...refSubcellPointProperties,
1037 typename paramPointValueType, class ...paramPointProperties>
1038 static void
1039 mapToReferenceSubcell( Kokkos::DynRankView<refSubcellPointValueType,refSubcellPointProperties...> refSubcellPoints,
1040 const Kokkos::DynRankView<paramPointValueType,paramPointProperties...> paramPoints,
1041 const ordinal_type subcellDim,
1042 const ordinal_type subcellOrd,
1043 const shards::CellTopology parentCell );
1044
1045
1046 //============================================================================================//
1047 // //
1048 // Physical-to-reference frame mapping and its inverse //
1049 // //
1050 //============================================================================================//
1051
1052
1096 template<typename refPointValueType, class ...refPointProperties,
1097 typename physPointValueType, class ...physPointProperties,
1098 typename worksetCellValueType, class ...worksetCellProperties>
1099 static void
1100 mapToReferenceFrame( Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
1101 const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1102 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1103 const shards::CellTopology cellTopo );
1104
1131 template<typename refPointValueType, class ...refPointProperties,
1132 typename initGuessValueType, class ...initGuessProperties,
1133 typename physPointValueType, class ...physPointProperties,
1134 typename worksetCellValueType, class ...worksetCellProperties,
1135 typename HGradBasisPtrType>
1136 static void
1137 mapToReferenceFrameInitGuess( Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
1138 const Kokkos::DynRankView<initGuessValueType,initGuessProperties...> initGuess,
1139 const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1140 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1141 const HGradBasisPtrType basis );
1142
1143
1174 template<typename refPointValueType, class ...refPointProperties,
1175 typename initGuessValueType, class ...initGuessProperties,
1176 typename physPointValueType, class ...physPointProperties,
1177 typename worksetCellValueType, class ...worksetCellProperties>
1178 static void
1179 mapToReferenceFrameInitGuess( Kokkos::DynRankView<refPointValueType,refPointProperties...> refPoints,
1180 const Kokkos::DynRankView<initGuessValueType,initGuessProperties...> initGuess,
1181 const Kokkos::DynRankView<physPointValueType,physPointProperties...> physPoints,
1182 const Kokkos::DynRankView<worksetCellValueType,worksetCellProperties...> worksetCell,
1183 const shards::CellTopology cellTopo ) {
1184 auto basis = createHGradBasis<refPointValueType,refPointValueType>(cellTopo);
1186 initGuess,
1187 physPoints,
1188 worksetCell,
1189 basis);
1190 }
1191
1192
1193 //============================================================================================//
1194 // //
1195 // Control Volume Coordinates //
1196 // //
1197 //============================================================================================//
1198
1272 template<typename subcvCoordValueType, class ...subcvCoordProperties,
1273 typename cellCoordValueType, class ...cellCoordProperties>
1274 static void
1275 getSubcvCoords( Kokkos::DynRankView<subcvCoordValueType,subcvCoordProperties...> subcvCoords,
1276 const Kokkos::DynRankView<cellCoordValueType,cellCoordProperties...> cellCoords,
1277 const shards::CellTopology primaryCell );
1278
1279 //============================================================================================//
1280 // //
1281 // Inclusion tests //
1282 // //
1283 //============================================================================================//
1284
1294 template<typename pointValueType, class ...pointProperties>
1295 static bool
1296 checkPointInclusion( const Kokkos::DynRankView<pointValueType,pointProperties...> point,
1297 const shards::CellTopology cellTopo,
1298 const double thres = threshold() );
1299
1300 // template<class ArrayPoint>
1301 // static int checkPointsetInclusion(const ArrayPoint & points,
1302 // const shards::CellTopology & cellTopo,
1303 // const double & threshold = INTREPID2_THRESHOLD);
1304
1305
1306
1333 template<typename inCellValueType, class ...inCellProperties,
1334 typename pointValueType, class ...pointProperties>
1335 static void checkPointwiseInclusion( Kokkos::DynRankView<inCellValueType,inCellProperties...> inCell,
1336 const Kokkos::DynRankView<pointValueType,pointProperties...> points,
1337 const shards::CellTopology cellTopo,
1338 const double thres = threshold() );
1339
1361 template<typename inCellValueType, class ...inCellProperties,
1362 typename pointValueType, class ...pointProperties,
1363 typename cellWorksetValueType, class ...cellWorksetProperties>
1364 static void checkPointwiseInclusion( Kokkos::DynRankView<inCellValueType,inCellProperties...> inCell,
1365 const Kokkos::DynRankView<pointValueType,pointProperties...> points,
1366 const Kokkos::DynRankView<cellWorksetValueType,cellWorksetProperties...> cellWorkset,
1367 const shards::CellTopology cellTopo,
1368 const double thres = threshold() );
1369
1370
1371 // //============================================================================================//
1372 // // //
1373 // // Debug //
1374 // // //
1375 // //============================================================================================//
1376
1377
1378 // /** \brief Prints the reference space coordinates of the vertices of the specified subcell
1379 // \param subcellDim [in] - dimension of the subcell where points are mapped to
1380 // \param subcellOrd [in] - subcell ordinal
1381 // \param parentCell [in] - cell topology of the parent cell.
1382 // */
1383 // static void printSubcellVertices(const int subcellDim,
1384 // const int subcellOrd,
1385 // const shards::CellTopology & parentCell);
1386
1387
1388
1389 // /** \brief Prints the nodes of a subcell from a cell workset
1390
1391 // */
1392 // template<class ArrayCell>
1393 // static void printWorksetSubcell(const ArrayCell & cellWorkset,
1394 // const shards::CellTopology & parentCell,
1395 // const int& pCellOrd,
1396 // const int& subcellDim,
1397 // const int& subcellOrd,
1398 // const int& fieldWidth = 3);
1399 };
1400
1401 //============================================================================================//
1402 // //
1403 // Validation of input/output arguments for CellTools methods //
1404 // //
1405 //============================================================================================//
1406
1413 template<typename jacobianViewType,
1414 typename PointViewType,
1415 typename worksetCellViewType>
1416 static void
1417 CellTools_setJacobianArgs( const jacobianViewType jacobian,
1418 const PointViewType points,
1419 const worksetCellViewType worksetCell,
1420 const shards::CellTopology cellTopo );
1421
1426 template<typename jacobianInvViewType,
1427 typename jacobianViewType>
1428 static void
1429 CellTools_setJacobianInvArgs( const jacobianInvViewType jacobianInv,
1430 const jacobianViewType jacobian );
1431
1432
1437 template<typename jacobianDetViewType,
1438 typename jacobianViewType>
1439 static void
1440 CellTools_setJacobianDetArgs( const jacobianDetViewType jacobianDet,
1441 const jacobianViewType jacobian );
1442
1443
1450 template<typename physPointViewType,
1451 typename refPointViewType,
1452 typename worksetCellViewType>
1453 static void
1454 CellTools_mapToPhysicalFrameArgs( const physPointViewType physPoints,
1455 const refPointViewType refPoints,
1456 const worksetCellViewType worksetCell,
1457 const shards::CellTopology cellTopo );
1458
1459
1466 template<typename refPointViewType,
1467 typename physPointViewType,
1468 typename worksetCellViewType>
1469 static void
1470 CellTools_mapToReferenceFrameArgs( const refPointViewType refPoints,
1471 const physPointViewType physPoints,
1472 const worksetCellViewType worksetCell,
1473 const shards::CellTopology cellTopo );
1474
1475
1476
1484 template<typename refPointViewType,
1485 typename initGuessViewType,
1486 typename physPointViewType,
1487 typename worksetCellViewType>
1488 static void
1489 CellTools_mapToReferenceFrameInitGuess( const refPointViewType refPoints,
1490 const initGuessViewType initGuess,
1491 const physPointViewType physPoints,
1492 const worksetCellViewType worksetCell,
1493 const shards::CellTopology cellTopo );
1494
1495 // /** \brief Validates arguments to Intrepid2::CellTools::checkPointwiseInclusion
1496 // \param inCell [out] - rank-1 (P) array required
1497 // \param physPoints [in] - rank-2 (P,D) array required
1498 // \param cellWorkset [in] - rank-3 (C,N,D) array required
1499 // \param whichCell [in] - 0 <= whichCell < C required
1500 // \param cellTopo [in] - cell topology with a reference cell required
1501 // */
1502 // template<class ArrayIncl, class ArrayPoint, class ArrayCell>
1503 // static void
1504 // validateArguments_checkPointwiseInclusion(ArrayIncl & inCell,
1505 // const ArrayPoint & physPoints,
1506 // const ArrayCell & cellWorkset,
1507 // const int & whichCell,
1508 // const shards::CellTopology & cell);
1509
1510
1511
1512}
1513
1515
1518
1522
1524
1525// not yet converted ...
1527
1528// #include "Intrepid2_CellToolsDefDebug.hpp"
1529
1530
1531#endif
1532
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_PYR_C1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_QUAD_C1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_QUAD_C2_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TET_C1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TET_C2_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TET_COMP12_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TRI_C1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_TRI_C2_FEM class.
Header file for the Intrepid2::Basis_HGRAD_WEDGE_C1_FEM class.
Header file for the Intrepid2::Basis_HGRAD_WEDGE_C2_FEM class.
Header file for Intrepid2::RealSpaceTools class providing basic linear algebra functionality in 1D,...
Contains definitions of custom data types in Intrepid2.
Header function for Intrepid2::Util class and other utility functions.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Hexahedron cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Hexahedron cell.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Line cell.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Pyramid cell.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Quadrilateral cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Quadrilateral cell.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Tetrahedron cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Tetrahedron cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Tetrahedron cell.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Triangle cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Triangle cell.
Implementation of the default H(grad)-compatible FEM basis of degree 1 on Wedge cell.
Implementation of the default H(grad)-compatible FEM basis of degree 2 on Wedge cell.
A stateless class for operations on cell data. Provides methods for:
static void 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 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 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.