49 #ifndef __INTREPID2_CELLTOOLS_DEF_NODE_INFO_HPP__
50 #define __INTREPID2_CELLTOOLS_DEF_NODE_INFO_HPP__
53 #if defined (__clang__) && !defined (__INTEL_COMPILER)
54 #pragma clang system_header
65 template<
typename SpT>
72 refNodeData_.line = referenceNodeDataViewType(
"CellTools::ReferenceNodeData::line", 2, 3);
73 refNodeData_.line_3 = referenceNodeDataViewType(
"CellTools::ReferenceNodeData::line_3", 3, 3);
75 refNodeData_.triangle = referenceNodeDataViewType(
"CellTools::ReferenceNodeData::triangle", 3, 3);
76 refNodeData_.triangle_4 = referenceNodeDataViewType(
"CellTools::ReferenceNodeData::triangle_4", 4, 3);
77 refNodeData_.triangle_6 = referenceNodeDataViewType(
"CellTools::ReferenceNodeData::triangle_6", 6, 3);
79 refNodeData_.quadrilateral = referenceNodeDataViewType(
"CellTools::ReferenceNodeData::quad", 4, 3);
80 refNodeData_.quadrilateral_8 = referenceNodeDataViewType(
"CellTools::ReferenceNodeData::quad_8", 8, 3);
81 refNodeData_.quadrilateral_9 = referenceNodeDataViewType(
"CellTools::ReferenceNodeData::quad_9", 9, 3);
83 refNodeData_.tetrahedron = referenceNodeDataViewType(
"CellTools::ReferenceNodeData::tet", 4, 3);
84 refNodeData_.tetrahedron_8 = referenceNodeDataViewType(
"CellTools::ReferenceNodeData::tet_8", 8, 3);
85 refNodeData_.tetrahedron_10 = referenceNodeDataViewType(
"CellTools::ReferenceNodeData::tet_10", 10, 3);
86 refNodeData_.tetrahedron_11 = referenceNodeDataViewType(
"CellTools::ReferenceNodeData::tet_11", 11, 3);
88 refNodeData_.hexahedron = referenceNodeDataViewType(
"CellTools::ReferenceNodeData::hex", 8, 3);
89 refNodeData_.hexahedron_20 = referenceNodeDataViewType(
"CellTools::ReferenceNodeData::hex_20", 20, 3);
90 refNodeData_.hexahedron_27 = referenceNodeDataViewType(
"CellTools::ReferenceNodeData::hex_27", 27, 3);
92 refNodeData_.pyramid = referenceNodeDataViewType(
"CellTools::ReferenceNodeData::pyr", 5, 3);
93 refNodeData_.pyramid_13 = referenceNodeDataViewType(
"CellTools::ReferenceNodeData::pyr_13", 13, 3);
94 refNodeData_.pyramid_14 = referenceNodeDataViewType(
"CellTools::ReferenceNodeData::pyr_14", 14, 3);
96 refNodeData_.wedge = referenceNodeDataViewType(
"CellTools::ReferenceNodeData::wedge", 6, 3);
97 refNodeData_.wedge_15 = referenceNodeDataViewType(
"CellTools::ReferenceNodeData::wedge_15", 15, 3);
98 refNodeData_.wedge_18 = referenceNodeDataViewType(
"CellTools::ReferenceNodeData::wedge_18", 18, 3);
103 Kokkos::deep_copy(refNodeData_.line, referenceNodeDataViewHostType(&refNodeDataStatic_.line[0][0], 2, 3));
104 Kokkos::deep_copy(refNodeData_.line_3, referenceNodeDataViewHostType(&refNodeDataStatic_.line_3[0][0], 3, 3));
106 Kokkos::deep_copy(refNodeData_.triangle, referenceNodeDataViewHostType(&refNodeDataStatic_.triangle[0][0], 3, 3));
107 Kokkos::deep_copy(refNodeData_.triangle_4, referenceNodeDataViewHostType(&refNodeDataStatic_.triangle_4[0][0], 4, 3));
108 Kokkos::deep_copy(refNodeData_.triangle_6, referenceNodeDataViewHostType(&refNodeDataStatic_.triangle_6[0][0], 6, 3));
110 Kokkos::deep_copy(refNodeData_.quadrilateral, referenceNodeDataViewHostType(&refNodeDataStatic_.quadrilateral[0][0], 4, 3));
111 Kokkos::deep_copy(refNodeData_.quadrilateral_8, referenceNodeDataViewHostType(&refNodeDataStatic_.quadrilateral_8[0][0], 8, 3));
112 Kokkos::deep_copy(refNodeData_.quadrilateral_9, referenceNodeDataViewHostType(&refNodeDataStatic_.quadrilateral_9[0][0], 9, 3));
114 Kokkos::deep_copy(refNodeData_.tetrahedron, referenceNodeDataViewHostType(&refNodeDataStatic_.tetrahedron[0][0], 4, 3));
115 Kokkos::deep_copy(refNodeData_.tetrahedron_8, referenceNodeDataViewHostType(&refNodeDataStatic_.tetrahedron_8[0][0], 8, 3));
116 Kokkos::deep_copy(refNodeData_.tetrahedron_10, referenceNodeDataViewHostType(&refNodeDataStatic_.tetrahedron_10[0][0], 10, 3));
117 Kokkos::deep_copy(refNodeData_.tetrahedron_11, referenceNodeDataViewHostType(&refNodeDataStatic_.tetrahedron_11[0][0], 11, 3));
119 Kokkos::deep_copy(refNodeData_.hexahedron, referenceNodeDataViewHostType(&refNodeDataStatic_.hexahedron[0][0], 8, 3));
120 Kokkos::deep_copy(refNodeData_.hexahedron_20, referenceNodeDataViewHostType(&refNodeDataStatic_.hexahedron_20[0][0], 20, 3));
121 Kokkos::deep_copy(refNodeData_.hexahedron_27, referenceNodeDataViewHostType(&refNodeDataStatic_.hexahedron_27[0][0], 27, 3));
123 Kokkos::deep_copy(refNodeData_.pyramid, referenceNodeDataViewHostType(&refNodeDataStatic_.pyramid[0][0], 5, 3));
124 Kokkos::deep_copy(refNodeData_.pyramid_13, referenceNodeDataViewHostType(&refNodeDataStatic_.pyramid_13[0][0], 13, 3));
125 Kokkos::deep_copy(refNodeData_.pyramid_14, referenceNodeDataViewHostType(&refNodeDataStatic_.pyramid_14[0][0], 14, 3));
127 Kokkos::deep_copy(refNodeData_.wedge, referenceNodeDataViewHostType(&refNodeDataStatic_.wedge[0][0], 6, 3));
128 Kokkos::deep_copy(refNodeData_.wedge_15, referenceNodeDataViewHostType(&refNodeDataStatic_.wedge_15[0][0], 15, 3));
129 Kokkos::deep_copy(refNodeData_.wedge_18, referenceNodeDataViewHostType(&refNodeDataStatic_.wedge_18[0][0], 18, 3));
132 Kokkos::push_finalize_hook( [=] {
133 refNodeData_.line = referenceNodeDataViewType();
134 refNodeData_.line_3 = referenceNodeDataViewType();
136 refNodeData_.triangle = referenceNodeDataViewType();
137 refNodeData_.triangle_4 = referenceNodeDataViewType();
138 refNodeData_.triangle_6 = referenceNodeDataViewType();
140 refNodeData_.quadrilateral = referenceNodeDataViewType();
141 refNodeData_.quadrilateral_8 = referenceNodeDataViewType();
142 refNodeData_.quadrilateral_9 = referenceNodeDataViewType();
144 refNodeData_.tetrahedron = referenceNodeDataViewType();
145 refNodeData_.tetrahedron_8 = referenceNodeDataViewType();
146 refNodeData_.tetrahedron_10 = referenceNodeDataViewType();
147 refNodeData_.tetrahedron_11 = referenceNodeDataViewType();
149 refNodeData_.hexahedron = referenceNodeDataViewType();
150 refNodeData_.hexahedron_20 = referenceNodeDataViewType();
151 refNodeData_.hexahedron_27 = referenceNodeDataViewType();
153 refNodeData_.pyramid = referenceNodeDataViewType();
154 refNodeData_.pyramid_13 = referenceNodeDataViewType();
155 refNodeData_.pyramid_14 = referenceNodeDataViewType();
157 refNodeData_.wedge = referenceNodeDataViewType();
158 refNodeData_.wedge_15 = referenceNodeDataViewType();
159 refNodeData_.wedge_18 = referenceNodeDataViewType();
162 isReferenceNodeDataSet_ =
true;
165 template<
typename SpT>
166 template<
typename cellCenterValueType,
class ...cellCenterProperties,
167 typename cellVertexValueType,
class ...cellVertexProperties>
171 Kokkos::DynRankView<cellVertexValueType,cellVertexProperties...> cellVertex,
172 const shards::CellTopology cell ) {
173 #ifdef HAVE_INTREPID2_DEBUG
174 INTREPID2_TEST_FOR_EXCEPTION( !hasReferenceCell(cell), std::invalid_argument,
175 ">>> ERROR (Intrepid2::CellTools::getReferenceCellCenter): the specified cell topology does not have a reference cell." );
177 INTREPID2_TEST_FOR_EXCEPTION( cellCenter.rank() != 1, std::invalid_argument,
178 ">>> ERROR (Intrepid2::CellTools::getReferenceCellCenter): cellCenter must have rank 1." );
180 INTREPID2_TEST_FOR_EXCEPTION( cellCenter.extent(0) < cell.getDimension(), std::invalid_argument,
181 ">>> ERROR (Intrepid2::CellTools::getReferenceCellCenter): cellCenter must have dimension bigger than Parameters::MaxDimension." );
183 INTREPID2_TEST_FOR_EXCEPTION( cellVertex.rank() != 1, std::invalid_argument,
184 ">>> ERROR (Intrepid2::CellTools::getReferenceCellCenter): cellVertex must have rank 1." );
186 INTREPID2_TEST_FOR_EXCEPTION( cellVertex.extent(0) < cell.getDimension(), std::invalid_argument,
187 ">>> ERROR (Intrepid2::CellTools::getReferenceCellCenter): cellVertex must have dimension bigger than Parameters::MaxDimension." );
189 const ordinal_type vertexCount = cell.getVertexCount();
190 const ordinal_type dim = cell.getDimension();
192 for (ordinal_type i=0;i<dim;++i) {
194 for (ordinal_type vertOrd=0;vertOrd<vertexCount;++vertOrd) {
195 getReferenceVertex(cellVertex, cell, vertOrd);
196 cellCenter(i) += cellVertex(i);
198 cellCenter(i) /= vertexCount;
203 template<
typename SpT>
204 template<
typename cellVertexValueType,
class ...cellVertexProperties>
208 const shards::CellTopology cell,
209 const ordinal_type vertexOrd ) {
210 #ifdef HAVE_INTREPID2_DEBUG
211 INTREPID2_TEST_FOR_EXCEPTION( !hasReferenceCell(cell), std::invalid_argument,
212 ">>> ERROR (Intrepid2::CellTools::getReferenceVertex): the specified cell topology does not have a reference cell." );
214 INTREPID2_TEST_FOR_EXCEPTION( (vertexOrd < 0) || vertexOrd > static_cast<ordinal_type>(cell.getVertexCount()), std::invalid_argument,
215 ">>> ERROR (Intrepid2::CellTools::getReferenceVertex): invalid node ordinal for the specified cell topology." );
217 INTREPID2_TEST_FOR_EXCEPTION( cellVertex.rank() != 1, std::invalid_argument,
218 ">>> ERROR (Intrepid2::CellTools::getReferenceNode): cellNodes must have rank 1." );
220 INTREPID2_TEST_FOR_EXCEPTION( cellVertex.extent(0) < cell.getDimension(), std::invalid_argument,
221 ">>> ERROR (Intrepid2::CellTools::getReferenceNode): cellNodes must have dimension bigger than Parameters::MaxDimension." );
223 getReferenceNode(cellVertex,
229 template<
typename SpT>
230 template<
typename subcellVertexValueType,
class ...subcellVertexProperties>
234 const ordinal_type subcellDim,
235 const ordinal_type subcellOrd,
236 const shards::CellTopology parentCell ) {
237 #ifdef HAVE_INTREPID2_DEBUG
238 INTREPID2_TEST_FOR_EXCEPTION( !hasReferenceCell(parentCell), std::invalid_argument,
239 ">>> ERROR (Intrepid2::CellTools::getReferenceSubcellVertices): the specified cell topology does not have a reference cell." );
241 INTREPID2_TEST_FOR_EXCEPTION( subcellDim > static_cast<ordinal_type>(parentCell.getDimension()), std::invalid_argument,
242 ">>> ERROR (Intrepid2::CellTools::getReferenceSubcellVertices): subcell dimension cannot exceed cell dimension." );
244 INTREPID2_TEST_FOR_EXCEPTION( subcellOrd >= static_cast<ordinal_type>(parentCell.getSubcellCount(subcellDim)), std::invalid_argument,
245 ">>> ERROR (Intrepid2::CellTools::getReferenceSubcellVertices): subcell ordinal cannot exceed subcell count." );
248 INTREPID2_TEST_FOR_EXCEPTION( subcellVertices.rank() != 2, std::invalid_argument,
249 ">>> ERROR (Intrepid2::CellTools::getReferenceSubcellVertices): subcellVertieces must have rank 2." );
252 INTREPID2_TEST_FOR_EXCEPTION( subcellVertices.extent(0) != parentCell.getNodeCount(subcellDim, subcellOrd), std::invalid_argument,
253 ">>> ERROR (Intrepid2::CellTools::getReferenceSubcellVertices): subcellVertieces dimension(0) must match to parent cell vertex count." );
255 INTREPID2_TEST_FOR_EXCEPTION( subcellVertices.extent(1) != parentCell.getDimension(), std::invalid_argument,
256 ">>> ERROR (Intrepid2::CellTools::getReferenceSubcellVertices): subcellVertieces dimension(1) must match to parent cell dimension." );
258 getReferenceSubcellNodes(subcellVertices,
265 template<
typename SpT>
266 template<
typename cellNodeValueType,
class ...cellNodeProperties>
270 const shards::CellTopology cell,
271 const ordinal_type nodeOrd ) {
272 #ifdef HAVE_INTREPID2_DEBUG
273 INTREPID2_TEST_FOR_EXCEPTION( !hasReferenceCell(cell), std::invalid_argument,
274 ">>> ERROR (Intrepid2::CellTools::getReferenceNode): the specified cell topology does not have a reference cell." );
276 INTREPID2_TEST_FOR_EXCEPTION( nodeOrd >= static_cast<ordinal_type>(cell.getNodeCount()), std::invalid_argument,
277 ">>> ERROR (Intrepid2::CellTools::getReferenceNode): invalid node ordinal for the specified cell topology." );
279 INTREPID2_TEST_FOR_EXCEPTION( cellNode.rank() != 1, std::invalid_argument,
280 ">>> ERROR (Intrepid2::CellTools::getReferenceNode): cellNode must have rank 1." );
282 INTREPID2_TEST_FOR_EXCEPTION( cellNode.extent(0) < cell.getDimension(), std::invalid_argument,
283 ">>> ERROR (Intrepid2::CellTools::getReferenceNode): cellNode must have dimension bigger than Parameters::MaxDimension." );
286 if (!isReferenceNodeDataSet_)
287 setReferenceNodeData();
289 referenceNodeDataViewType ref;
291 switch (cell.getKey() ) {
292 case shards::Line<2>::key:
293 case shards::ShellLine<2>::key:
294 case shards::Beam<2>::key: ref = refNodeData_.line;
break;
295 case shards::Line<3>::key:
296 case shards::ShellLine<3>::key:
297 case shards::Beam<3>::key: ref = refNodeData_.line_3;
break;
299 case shards::Triangle<3>::key:
300 case shards::ShellTriangle<3>::key: ref = refNodeData_.triangle;
break;
301 case shards::Triangle<4>::key: ref = refNodeData_.triangle_4;
break;
302 case shards::Triangle<6>::key:
303 case shards::ShellTriangle<6>::key: ref = refNodeData_.triangle_6;
break;
305 case shards::Quadrilateral<4>::key:
306 case shards::ShellQuadrilateral<4>::key: ref = refNodeData_.quadrilateral;
break;
307 case shards::Quadrilateral<8>::key:
308 case shards::ShellQuadrilateral<8>::key: ref = refNodeData_.quadrilateral_8;
break;
309 case shards::Quadrilateral<9>::key:
310 case shards::ShellQuadrilateral<9>::key: ref = refNodeData_.quadrilateral_9;
break;
312 case shards::Tetrahedron<4>::key: ref = refNodeData_.tetrahedron;
break;
313 case shards::Tetrahedron<8>::key: ref = refNodeData_.tetrahedron_8;
break;
314 case shards::Tetrahedron<10>::key: ref = refNodeData_.tetrahedron_10;
break;
315 case shards::Tetrahedron<11>::key: ref = refNodeData_.tetrahedron_11;
break;
317 case shards::Hexahedron<8>::key: ref = refNodeData_.hexahedron;
break;
318 case shards::Hexahedron<20>::key: ref = refNodeData_.hexahedron_20;
break;
319 case shards::Hexahedron<27>::key: ref = refNodeData_.hexahedron_27;
break;
321 case shards::Pyramid<5>::key: ref = refNodeData_.pyramid;
break;
322 case shards::Pyramid<13>::key: ref = refNodeData_.pyramid_13;
break;
323 case shards::Pyramid<14>::key: ref = refNodeData_.pyramid_14;
break;
325 case shards::Wedge<6>::key: ref = refNodeData_.wedge;
break;
326 case shards::Wedge<15>::key: ref = refNodeData_.wedge_15;
break;
327 case shards::Wedge<18>::key: ref = refNodeData_.wedge_18;
break;
330 INTREPID2_TEST_FOR_EXCEPTION(
true, std::invalid_argument,
331 ">>> ERROR (Intrepid2::CellTools::getReferenceNode): invalid cell topology.");
339 const ordinal_type dim = cell.getDimension();
341 for (ordinal_type i=0;i<dim;++i)
342 cellNode(i) = ref(nodeOrd, i);
345 template<
typename SpT>
346 template<
typename subcellNodeValueType,
class ...subcellNodeProperties>
350 const ordinal_type subcellDim,
351 const ordinal_type subcellOrd,
352 const shards::CellTopology parentCell ) {
353 #ifdef HAVE_INTREPID2_DEBUG
354 INTREPID2_TEST_FOR_EXCEPTION( !hasReferenceCell(parentCell), std::invalid_argument,
355 ">>> ERROR (Intrepid2::CellTools::getReferenceSubcellNodes): the specified cell topology does not have a reference cell.");
357 INTREPID2_TEST_FOR_EXCEPTION( subcellDim > static_cast<ordinal_type>(parentCell.getDimension()), std::invalid_argument,
358 ">>> ERROR (Intrepid2::CellTools::getReferenceSubcellNodes): subcell dimension out of range.");
360 INTREPID2_TEST_FOR_EXCEPTION( subcellOrd >= static_cast<ordinal_type>(parentCell.getSubcellCount(subcellDim)), std::invalid_argument,
361 ">>> ERROR (Intrepid2::CellTools::getReferenceSubcellNodes): subcell ordinal out of range.");
364 INTREPID2_TEST_FOR_EXCEPTION( subcellNodes.rank() != 2, std::invalid_argument,
365 ">>> ERROR (Intrepid2::CellTools::getReferenceSubcellNodes): subcellNodes must have rank 2.");
367 INTREPID2_TEST_FOR_EXCEPTION( subcellNodes.extent(0) != parentCell.getNodeCount(subcellDim, subcellOrd), std::invalid_argument,
368 ">>> ERROR (Intrepid2::CellTools::getReferenceSubcellNodes): subcellNodes dimension (0) must match to node count in cell.");
370 INTREPID2_TEST_FOR_EXCEPTION( subcellNodes.extent(1) != parentCell.getDimension(), std::invalid_argument,
371 ">>> ERROR (Intrepid2::CellTools::getReferenceSubcellNodes): subcellNodes dimension (1) must match to cell dimension.");
375 const auto subcNodeCount = parentCell.getNodeCount(subcellDim, subcellOrd);
378 for (size_type subcNodeOrd=0;subcNodeOrd<subcNodeCount;++subcNodeOrd) {
380 const auto cellNodeOrd = parentCell.getNodeMap(subcellDim, subcellOrd, subcNodeOrd);
382 auto dst = Kokkos::subdynrankview(subcellNodes, subcNodeOrd, Kokkos::ALL());
383 getReferenceNode(dst, parentCell, cellNodeOrd);
387 template<
typename SpT>
388 template<
typename refEdgeTangentValueType,
class ...refEdgeTangentProperties>
392 const ordinal_type edgeOrd,
393 const shards::CellTopology parentCell ) {
394 #ifdef HAVE_INTREPID2_DEBUG
395 INTREPID2_TEST_FOR_EXCEPTION( parentCell.getDimension() != 2 &&
396 parentCell.getDimension() != 3, std::invalid_argument,
397 ">>> ERROR (Intrepid2::CellTools::getReferenceFaceTangents): two or three-dimensional parent cell required");
399 INTREPID2_TEST_FOR_EXCEPTION( refEdgeTangent.rank() != 1, std::invalid_argument,
400 ">>> ERROR (Intrepid2::CellTools::getReferenceFaceTangents): rank = 1 required for output arrays");
402 INTREPID2_TEST_FOR_EXCEPTION( refEdgeTangent.extent(0) != parentCell.getDimension(), std::invalid_argument,
403 ">>> ERROR (Intrepid2::CellTools::getReferenceFaceTangents): output array size is required to match space dimension");
405 INTREPID2_TEST_FOR_EXCEPTION( edgeOrd < 0 ||
406 edgeOrd >= static_cast<ordinal_type>(parentCell.getSubcellCount(1)), std::invalid_argument,
407 ">>> ERROR (Intrepid2::CellTools::getReferenceFaceTangents): edge ordinal out of bounds");
412 subcellParamViewType edgeMap;
413 getSubcellParametrization(edgeMap, 1, parentCell);
417 const ordinal_type dim = parentCell.getDimension();
418 for (ordinal_type i=0;i<dim;++i)
419 refEdgeTangent(i) = edgeMap(edgeOrd, i, 1);
423 template<
typename SpT>
424 template<
typename refFaceTanUValueType,
class ...refFaceTanUProperties,
425 typename refFaceTanVValueType,
class ...refFaceTanVProperties>
429 Kokkos::DynRankView<refFaceTanVValueType,refFaceTanVProperties...> refFaceTanV,
430 const ordinal_type faceOrd,
431 const shards::CellTopology parentCell ) {
432 #ifdef HAVE_INTREPID2_DEBUG
433 INTREPID2_TEST_FOR_EXCEPTION( parentCell.getDimension() != 3, std::invalid_argument,
434 ">>> ERROR (Intrepid2::CellTools::getReferenceFaceTangents): three-dimensional parent cell required");
436 INTREPID2_TEST_FOR_EXCEPTION( faceOrd < 0 || faceOrd >= static_cast<ordinal_type>(parentCell.getSubcellCount(2)), std::invalid_argument,
437 ">>> ERROR (Intrepid2::CellTools::getReferenceFaceTangents): face ordinal out of bounds");
439 INTREPID2_TEST_FOR_EXCEPTION( refFaceTanU.rank() != 1 || refFaceTanV.rank() != 1, std::invalid_argument,
440 ">>> ERROR (Intrepid2::CellTools::getReferenceFaceTangents): rank = 1 required for output arrays");
442 INTREPID2_TEST_FOR_EXCEPTION( refFaceTanU.extent(0) != parentCell.getDimension(), std::invalid_argument,
443 ">>> ERROR (Intrepid2::CellTools::getReferenceFaceTangents): dim0 (spatial dim) must match that of parent cell");
445 INTREPID2_TEST_FOR_EXCEPTION( refFaceTanV.extent(0) != parentCell.getDimension(), std::invalid_argument,
446 ">>> ERROR (Intrepid2::CellTools::getReferenceFaceTangents): dim0 (spatial dim) must match that of parent cell");
451 subcellParamViewType faceMap;
452 getSubcellParametrization(faceMap, 2, parentCell);
460 const ordinal_type dim = parentCell.getDimension();
461 for (ordinal_type i=0;i<dim;++i) {
462 refFaceTanU(i) = faceMap(faceOrd, i, 1);
463 refFaceTanV(i) = faceMap(faceOrd, i, 2);
467 template<
typename SpT>
468 template<
typename refSideNormalValueType,
class ...refSideNormalProperties>
472 const ordinal_type sideOrd,
473 const shards::CellTopology parentCell ) {
474 #ifdef HAVE_INTREPID2_DEBUG
475 INTREPID2_TEST_FOR_EXCEPTION( parentCell.getDimension() != 2 &&
476 parentCell.getDimension() != 3, std::invalid_argument,
477 ">>> ERROR (Intrepid2::CellTools::getReferenceSideNormal): two or three-dimensional parent cell required");
480 INTREPID2_TEST_FOR_EXCEPTION( sideOrd < 0 || sideOrd >= static_cast<ordinal_type>(parentCell.getSubcellCount(parentCell.getDimension() - 1)), std::invalid_argument,
481 ">>> ERROR (Intrepid2::CellTools::getReferenceSideNormal): side ordinal out of bounds");
484 const auto dim = parentCell.getDimension();
487 getReferenceEdgeTangent(refSideNormal, sideOrd, parentCell);
490 refSideNormalValueType tmp = refSideNormal(0);
491 refSideNormal(0) = refSideNormal(1);
492 refSideNormal(1) = -tmp;
495 getReferenceFaceNormal(refSideNormal, sideOrd, parentCell);
500 template<
typename SpT>
501 template<
typename refFaceNormalValueType,
class ...refFaceNormalProperties>
505 const ordinal_type faceOrd,
506 const shards::CellTopology parentCell ) {
507 #ifdef HAVE_INTREPID2_DEBUG
508 INTREPID2_TEST_FOR_EXCEPTION( parentCell.getDimension() != 3, std::invalid_argument,
509 ">>> ERROR (Intrepid2::CellTools::getReferenceFaceNormal): three-dimensional parent cell required");
511 INTREPID2_TEST_FOR_EXCEPTION( faceOrd < 0 || faceOrd >= static_cast<ordinal_type>(parentCell.getSubcellCount(2)), std::invalid_argument,
512 ">>> ERROR (Intrepid2::CellTools::getReferenceFaceNormal): face ordinal out of bounds");
514 INTREPID2_TEST_FOR_EXCEPTION( refFaceNormal.rank() != 1, std::invalid_argument,
515 ">>> ERROR (Intrepid2::CellTools::getReferenceFaceNormal): rank = 1 required for output array");
517 INTREPID2_TEST_FOR_EXCEPTION( refFaceNormal.extent(0) != parentCell.getDimension(), std::invalid_argument,
518 ">>> ERROR (Intrepid2::CellTools::getReferenceFaceNormal): dim0 (spatial dim) must match that of parent cell");
522 const auto dim = parentCell.getDimension();
523 auto vcprop = Kokkos::common_view_alloc_prop(refFaceNormal);
524 using common_value_type =
typename decltype(vcprop)::value_type;
525 Kokkos::DynRankView< common_value_type, SpT > refFaceTanU ( Kokkos::view_alloc(
"CellTools::getReferenceFaceNormal::refFaceTanU", vcprop), dim );
526 Kokkos::DynRankView< common_value_type, SpT > refFaceTanV ( Kokkos::view_alloc(
"CellTools::getReferenceFaceNormal::refFaceTanV", vcprop), dim );
527 getReferenceFaceTangents(refFaceTanU, refFaceTanV, faceOrd, parentCell);
532 template<
typename SpT>
533 template<
typename edgeTangentValueType,
class ...edgeTangentProperties,
534 typename worksetJacobianValueType,
class ...worksetJacobianProperties>
538 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
539 const ordinal_type worksetEdgeOrd,
540 const shards::CellTopology parentCell ) {
541 #ifdef HAVE_INTREPID2_DEBUG
542 INTREPID2_TEST_FOR_EXCEPTION( parentCell.getDimension() != 3 &&
543 parentCell.getDimension() != 2, std::invalid_argument,
544 ">>> ERROR (Intrepid2::CellTools::getPhysicalEdgeTangents): 2D or 3D parent cell required." );
547 INTREPID2_TEST_FOR_EXCEPTION( edgeTangents.rank() != 3, std::invalid_argument,
548 ">>> ERROR (Intrepid2::CellTools::getPhysicalEdgeTangents): edgeTangents requires rank 3." );
549 INTREPID2_TEST_FOR_EXCEPTION( edgeTangents.extent(2) != 2 &&
550 edgeTangents.extent(2) != 3, std::invalid_argument,
551 ">>> ERROR (Intrepid2::CellTools::getPhysicalEdgeTangents): edgeTangents dimension(2) must be 2 or 3." );
554 INTREPID2_TEST_FOR_EXCEPTION( worksetJacobians.rank() != 4, std::invalid_argument,
555 ">>> ERROR (Intrepid2::CellTools::getPhysicalEdgeTangents): worksetJacobians requires rank 4." );
556 INTREPID2_TEST_FOR_EXCEPTION( worksetJacobians.extent(2) != 2 &&
557 worksetJacobians.extent(2) != 3, std::invalid_argument,
558 ">>> ERROR (Intrepid2::CellTools::getPhysicalEdgeTangents): worksetJacobians dimension(2) must be 2 or 3." );
559 INTREPID2_TEST_FOR_EXCEPTION( worksetJacobians.extent(2) != worksetJacobians.extent(3), std::invalid_argument,
560 ">>> ERROR (Intrepid2::CellTools::getPhysicalEdgeTangents): worksetJacobians dimension(2) and (3) must match each other." );
563 for (
auto i=0;i<3;++i) {
564 INTREPID2_TEST_FOR_EXCEPTION( edgeTangents.extent(i) != worksetJacobians.extent(i), std::invalid_argument,
565 ">>> ERROR (Intrepid2::CellTools::getPhysicalEdgeTangents): edgeTangents dimension (i) does not match to worksetJacobians dimension(i)." );
570 const auto dim = parentCell.getDimension();
571 auto vcprop = Kokkos::common_view_alloc_prop(edgeTangents);
572 using common_value_type =
typename decltype(vcprop)::value_type;
573 Kokkos::DynRankView< common_value_type, SpT > refEdgeTan ( Kokkos::view_alloc(
"CellTools::getPhysicalEdgeTangents::refEdgeTan", vcprop), dim);
574 getReferenceEdgeTangent(refEdgeTan, worksetEdgeOrd, parentCell);
580 template<
typename SpT>
581 template<
typename faceTanUValueType,
class ...faceTanUProperties,
582 typename faceTanVValueType,
class ...faceTanVProperties,
583 typename worksetJacobianValueType,
class ...worksetJacobianProperties>
587 Kokkos::DynRankView<faceTanVValueType,faceTanVProperties...> faceTanV,
588 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
589 const ordinal_type worksetFaceOrd,
590 const shards::CellTopology parentCell ) {
591 #ifdef HAVE_INTREPID2_DEBUG
592 INTREPID2_TEST_FOR_EXCEPTION( parentCell.getDimension() != 3, std::invalid_argument,
593 ">>> ERROR (Intrepid2::CellTools::getPhysicalFaceTangents): three-dimensional parent cell required");
596 INTREPID2_TEST_FOR_EXCEPTION( faceTanU.rank() != 3 ||
597 faceTanV.rank() != 3, std::invalid_argument,
598 ">>> ERROR (Intrepid2::CellTools::getPhysicalFaceTangents): faceTan U,V must have rank 3." );
600 INTREPID2_TEST_FOR_EXCEPTION( faceTanU.extent(2) != 3 ||
601 faceTanV.extent(2) != 3, std::invalid_argument,
602 ">>> ERROR (Intrepid2::CellTools::getPhysicalFaceTangents): faceTan U,V dimension (2) must be 3." );
604 for (
auto i=0;i<3;++i) {
605 INTREPID2_TEST_FOR_EXCEPTION( faceTanU.extent(i) != faceTanV.extent(i), std::invalid_argument,
606 ">>> ERROR (Intrepid2::CellTools::getPhysicalFaceTangents): faceTan U,V dimension (i) must match each other." );
610 INTREPID2_TEST_FOR_EXCEPTION( worksetJacobians.rank() != 4, std::invalid_argument,
611 ">>> ERROR (Intrepid2::CellTools::getPhysicalFaceTangents): worksetJacobians must have rank 4." );
613 INTREPID2_TEST_FOR_EXCEPTION( worksetJacobians.extent(2) != 3, std::invalid_argument,
614 ">>> ERROR (Intrepid2::CellTools::getPhysicalFaceTangents): worksetJacobians dimension(2) must be 3." );
616 INTREPID2_TEST_FOR_EXCEPTION( worksetJacobians.extent(2) != worksetJacobians.extent(3), std::invalid_argument,
617 ">>> ERROR (Intrepid2::CellTools::getPhysicalFaceTangents): worksetJacobians dimension(2) and dimension(3) must match." );
620 for (
auto i=0;i<3;++i) {
621 INTREPID2_TEST_FOR_EXCEPTION( faceTanU.extent(i) != worksetJacobians.extent(i), std::invalid_argument,
622 ">>> ERROR (Intrepid2::CellTools::getPhysicalFaceTangents): worksetJacobians dimension(i) and faceTan dimension (i) must match." );
627 const auto dim = parentCell.getDimension();
629 auto vcpropU = Kokkos::common_view_alloc_prop(faceTanU);
630 using common_value_typeU =
typename decltype(vcpropU)::value_type;
631 Kokkos::DynRankView< common_value_typeU, SpT > refFaceTanU ( Kokkos::view_alloc(
"CellTools::getPhysicalFaceTangents::refFaceTanU", vcpropU), dim);
633 auto vcpropV = Kokkos::common_view_alloc_prop(faceTanV);
634 using common_value_typeV =
typename decltype(vcpropV)::value_type;
635 Kokkos::DynRankView< common_value_typeV, SpT > refFaceTanV ( Kokkos::view_alloc(
"CellTools::getPhysicalFaceTangents::refFaceTanV", vcpropV), dim);
637 getReferenceFaceTangents(refFaceTanU, refFaceTanV, worksetFaceOrd, parentCell);
644 template<
typename SpT>
645 template<
typename sideNormalValueType,
class ...sideNormalProperties,
646 typename worksetJacobianValueType,
class ...worksetJacobianProperties>
650 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
651 const ordinal_type worksetSideOrd,
652 const shards::CellTopology parentCell ) {
653 #ifdef HAVE_INTREPID2_DEBUG
654 INTREPID2_TEST_FOR_EXCEPTION( parentCell.getDimension() != 2 &&
655 parentCell.getDimension() != 3, std::invalid_argument,
656 ">>> ERROR (Intrepid2::CellTools::getPhysicalSideNormals): two or three-dimensional parent cell required");
659 INTREPID2_TEST_FOR_EXCEPTION( worksetSideOrd < 0 ||
660 worksetSideOrd >= static_cast<ordinal_type>(parentCell.getSubcellCount(parentCell.getDimension() - 1)), std::invalid_argument,
661 ">>> ERROR (Intrepid2::CellTools::getPhysicalSideNormals): side ordinal out of bounds");
663 const auto dim = parentCell.getDimension();
667 auto vcprop = Kokkos::common_view_alloc_prop(sideNormals);
668 using common_value_type =
typename decltype(vcprop)::value_type;
669 Kokkos::DynRankView< common_value_type, SpT > edgeTangents ( Kokkos::view_alloc(
"CellTools::getPhysicalSideNormals::edgeTan", vcprop),
670 sideNormals.extent(0),
671 sideNormals.extent(1),
672 sideNormals.extent(2));
673 getPhysicalEdgeTangents(edgeTangents, worksetJacobians, worksetSideOrd, parentCell);
675 Kokkos::DynRankView< common_value_type, SpT > rotation ( Kokkos::view_alloc(
"CellTools::getPhysicalSideNormals::rotation", vcprop), dim, dim);
676 rotation(0,0) = 0; rotation(0,1) = 1;
677 rotation(1,0) = -1; rotation(1,1) = 0;
681 getPhysicalFaceNormals(sideNormals, worksetJacobians, worksetSideOrd, parentCell);
686 template<
typename SpT>
687 template<
typename faceNormalValueType,
class ...faceNormalProperties,
688 typename worksetJacobianValueType,
class ...worksetJacobianProperties>
692 const Kokkos::DynRankView<worksetJacobianValueType,worksetJacobianProperties...> worksetJacobians,
693 const ordinal_type worksetFaceOrd,
694 const shards::CellTopology parentCell ) {
695 #ifdef HAVE_INTREPID2_DEBUG
696 INTREPID2_TEST_FOR_EXCEPTION( parentCell.getDimension() != 3, std::invalid_argument,
697 ">>> ERROR (Intrepid2::CellTools::getPhysicalFaceNormals): three-dimensional parent cell required." );
700 INTREPID2_TEST_FOR_EXCEPTION( faceNormals.rank() != 3, std::invalid_argument,
701 ">>> ERROR (Intrepid2::CellTools::getPhysicalFaceNormals): faceNormals must have a rank 3." );
702 INTREPID2_TEST_FOR_EXCEPTION( faceNormals.extent(2) != 3, std::invalid_argument,
703 ">>> ERROR (Intrepid2::CellTools::getPhysicalFaceNormals): faceNormals dimension (2) must be 3." );
706 INTREPID2_TEST_FOR_EXCEPTION( worksetJacobians.rank() != 4, std::invalid_argument,
707 ">>> ERROR (Intrepid2::CellTools::getPhysicalFaceNormals): worksetJacobians must have a rank 4." );
708 INTREPID2_TEST_FOR_EXCEPTION( worksetJacobians.extent(2) != 3, std::invalid_argument,
709 ">>> ERROR (Intrepid2::CellTools::getPhysicalFaceNormals): worksetJacobians dimension (2) must be 3." );
710 INTREPID2_TEST_FOR_EXCEPTION( worksetJacobians.extent(2) != worksetJacobians.extent(3), std::invalid_argument,
711 ">>> ERROR (Intrepid2::CellTools::getPhysicalFaceNormals): worksetJacobians dimension (2) must match to dimension (3)." );
714 for (
auto i=0;i<3;++i) {
715 INTREPID2_TEST_FOR_EXCEPTION( faceNormals.extent(i) != worksetJacobians.extent(i), std::invalid_argument,
716 ">>> ERROR (Intrepid2::CellTools::getPhysicalFaceNormals): faceNormals dimension (i) must match to worksetJacobians dimension (i)." );
722 const auto worksetSize = worksetJacobians.extent(0);
723 const auto facePtCount = worksetJacobians.extent(1);
724 const auto dim = parentCell.getDimension();
726 auto vcprop = Kokkos::common_view_alloc_prop(faceNormals);
727 using common_value_type =
typename decltype(vcprop)::value_type;
728 Kokkos::DynRankView< common_value_type, SpT > faceTanU ( Kokkos::view_alloc(
"CellTools::getPhysicalFaceNormals::faceTanU", vcprop), worksetSize, facePtCount, dim);
729 Kokkos::DynRankView< common_value_type, SpT > faceTanV ( Kokkos::view_alloc(
"CellTools::getPhysicalFaceNormals::faceTanV", vcprop), worksetSize, facePtCount, dim);
731 getPhysicalFaceTangents(faceTanU, faceTanV,
740 template<
typename SpT>
743 isReferenceNodeDataSet_ =
false;
745 template<
typename SpT>
750 template<
typename SpT>
753 refNodeDataStatic_ = {
756 {-1.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}
759 {-1.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 0.0, 0.0}
763 { 0.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0}
766 { 0.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0}, { 1/3, 1/3, 0.0}
769 { 0.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0},
770 { 0.5, 0.0, 0.0}, { 0.5, 0.5, 0.0}, { 0.0, 0.5, 0.0}
774 {-1.0,-1.0, 0.0}, { 1.0,-1.0, 0.0}, { 1.0, 1.0, 0.0}, {-1.0, 1.0, 0.0}
777 {-1.0,-1.0, 0.0}, { 1.0,-1.0, 0.0}, { 1.0, 1.0, 0.0}, {-1.0, 1.0, 0.0},
778 { 0.0,-1.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0}, {-1.0, 0.0, 0.0}
781 {-1.0,-1.0, 0.0}, { 1.0,-1.0, 0.0}, { 1.0, 1.0, 0.0}, {-1.0, 1.0, 0.0},
782 { 0.0,-1.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0}, {-1.0, 0.0, 0.0}, { 0.0, 0.0, 0.0}
786 { 0.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0}, { 0.0, 0.0, 1.0}
789 { 0.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0}, { 0.0, 0.0, 1.0},
790 { 1/3, 0.0, 1/3}, { 1/3, 1/3, 1/3}, { 1/3, 1/3, 0.0}, { 0.0, 1/3, 1/3}
793 { 0.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0}, { 0.0, 0.0, 1.0},
794 { 0.5, 0.0, 0.0}, { 0.5, 0.5, 0.0}, { 0.0, 0.5, 0.0}, { 0.0, 0.0, 0.5}, { 0.5, 0.0, 0.5}, { 0.0, 0.5, 0.5}
797 { 0.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0}, { 0.0, 0.0, 1.0},
798 { 0.5, 0.0, 0.0}, { 0.5, 0.5, 0.0}, { 0.0, 0.5, 0.0}, { 0.0, 0.0, 0.5}, { 0.5, 0.0, 0.5}, { 0.0, 0.5, 0.5}
802 {-1.0,-1.0,-1.0}, { 1.0,-1.0,-1.0}, { 1.0, 1.0,-1.0}, {-1.0, 1.0,-1.0},
803 {-1.0,-1.0, 1.0}, { 1.0,-1.0, 1.0}, { 1.0, 1.0, 1.0}, {-1.0, 1.0, 1.0}
806 {-1.0,-1.0,-1.0}, { 1.0,-1.0,-1.0}, { 1.0, 1.0,-1.0}, {-1.0, 1.0,-1.0},
807 {-1.0,-1.0, 1.0}, { 1.0,-1.0, 1.0}, { 1.0, 1.0, 1.0}, {-1.0, 1.0, 1.0},
808 { 0.0,-1.0,-1.0}, { 1.0, 0.0,-1.0}, { 0.0, 1.0,-1.0}, {-1.0, 0.0,-1.0},
809 {-1.0,-1.0, 0.0}, { 1.0,-1.0, 0.0}, { 1.0, 1.0, 0.0}, {-1.0, 1.0, 0.0},
810 { 0.0,-1.0, 1.0}, { 1.0, 0.0, 1.0}, { 0.0, 1.0, 1.0}, {-1.0, 0.0, 1.0}
813 {-1.0,-1.0,-1.0}, { 1.0,-1.0,-1.0}, { 1.0, 1.0,-1.0}, {-1.0, 1.0,-1.0},
814 {-1.0,-1.0, 1.0}, { 1.0,-1.0, 1.0}, { 1.0, 1.0, 1.0}, {-1.0, 1.0, 1.0},
815 { 0.0,-1.0,-1.0}, { 1.0, 0.0,-1.0}, { 0.0, 1.0,-1.0}, {-1.0, 0.0,-1.0},
816 {-1.0,-1.0, 0.0}, { 1.0,-1.0, 0.0}, { 1.0, 1.0, 0.0}, {-1.0, 1.0, 0.0},
817 { 0.0,-1.0, 1.0}, { 1.0, 0.0, 1.0}, { 0.0, 1.0, 1.0}, {-1.0, 0.0, 1.0},
819 { 0.0, 0.0,-1.0}, { 0.0, 0.0, 1.0}, {-1.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}, {0.0,-1.0, 0.0}, {0.0, 1.0, 0.0}
823 {-1.0,-1.0, 0.0}, { 1.0,-1.0, 0.0}, { 1.0, 1.0, 0.0}, {-1.0, 1.0, 0.0}, { 0.0, 0.0, 1.0}
826 {-1.0,-1.0, 0.0}, { 1.0,-1.0, 0.0}, { 1.0, 1.0, 0.0}, {-1.0, 1.0, 0.0}, { 0.0, 0.0, 1.0},
827 { 0.0,-1.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0}, {-1.0, 0.0, 0.0},
828 {-0.5,-0.5, 0.5}, { 0.5,-0.5, 0.5}, { 0.5, 0.5, 0.5}, {-0.5, 0.5, 0.5}
831 {-1.0,-1.0, 0.0}, { 1.0,-1.0, 0.0}, { 1.0, 1.0, 0.0}, {-1.0, 1.0, 0.0}, { 0.0, 0.0, 1.0},
832 { 0.0,-1.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0}, {-1.0, 0.0, 0.0},
833 {-0.5,-0.5, 0.5}, { 0.5,-0.5, 0.5}, { 0.5, 0.5, 0.5}, {-0.5, 0.5, 0.5}, { 0.0, 0.0, 0.0}
837 { 0.0, 0.0,-1.0}, { 1.0, 0.0,-1.0}, { 0.0, 1.0,-1.0}, { 0.0, 0.0, 1.0}, { 1.0, 0.0, 1.0}, { 0.0, 1.0, 1.0}
840 { 0.0, 0.0,-1.0}, { 1.0, 0.0,-1.0}, { 0.0, 1.0,-1.0}, { 0.0, 0.0, 1.0}, { 1.0, 0.0, 1.0}, { 0.0, 1.0, 1.0},
841 { 0.5, 0.0,-1.0}, { 0.5, 0.5,-1.0}, { 0.0, 0.5,-1.0}, { 0.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0},
842 { 0.5, 0.0, 1.0}, { 0.5, 0.5, 1.0}, { 0.0, 0.5, 1.0}
845 { 0.0, 0.0,-1.0}, { 1.0, 0.0,-1.0}, { 0.0, 1.0,-1.0}, { 0.0, 0.0, 1.0}, { 1.0, 0.0, 1.0}, { 0.0, 1.0, 1.0},
846 { 0.5, 0.0,-1.0}, { 0.5, 0.5,-1.0}, { 0.0, 0.5,-1.0}, { 0.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0},
847 { 0.5, 0.0, 1.0}, { 0.5, 0.5, 1.0}, { 0.0, 0.5, 1.0},
848 { 0.5, 0.0, 0.0}, { 0.5, 0.5, 0.0}, { 0.0, 0.5, 0.0}