Intrepid2
Intrepid2_CellToolsDefNodeInfo.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 
43 
49 #ifndef __INTREPID2_CELLTOOLS_DEF_NODE_INFO_HPP__
50 #define __INTREPID2_CELLTOOLS_DEF_NODE_INFO_HPP__
51 
52 // disable clang warnings
53 #if defined (__clang__) && !defined (__INTEL_COMPILER)
54 #pragma clang system_header
55 #endif
56 
57 namespace Intrepid2 {
58 
59  //============================================================================================//
60  // //
61  // Reference nodes //
62  // //
63  //============================================================================================//
64 
65  template<typename SpT>
66  void
69 
70  {
71  // create memory on devices
72  refNodeData_.line = referenceNodeDataViewType("CellTools::ReferenceNodeData::line", 2, 3);
73  refNodeData_.line_3 = referenceNodeDataViewType("CellTools::ReferenceNodeData::line_3", 3, 3);
74 
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);
78 
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);
82 
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);
87 
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);
91 
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);
95 
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);
99  }
100 
101  {
102  // copy static data to devices
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));
105 
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));
109 
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));
113 
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));
118 
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));
122 
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));
126 
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));
130  }
131 
132  Kokkos::push_finalize_hook( [=] {
133  refNodeData_.line = referenceNodeDataViewType();
134  refNodeData_.line_3 = referenceNodeDataViewType();
135 
136  refNodeData_.triangle = referenceNodeDataViewType();
137  refNodeData_.triangle_4 = referenceNodeDataViewType();
138  refNodeData_.triangle_6 = referenceNodeDataViewType();
139 
140  refNodeData_.quadrilateral = referenceNodeDataViewType();
141  refNodeData_.quadrilateral_8 = referenceNodeDataViewType();
142  refNodeData_.quadrilateral_9 = referenceNodeDataViewType();
143 
144  refNodeData_.tetrahedron = referenceNodeDataViewType();
145  refNodeData_.tetrahedron_8 = referenceNodeDataViewType();
146  refNodeData_.tetrahedron_10 = referenceNodeDataViewType();
147  refNodeData_.tetrahedron_11 = referenceNodeDataViewType();
148 
149  refNodeData_.hexahedron = referenceNodeDataViewType();
150  refNodeData_.hexahedron_20 = referenceNodeDataViewType();
151  refNodeData_.hexahedron_27 = referenceNodeDataViewType();
152 
153  refNodeData_.pyramid = referenceNodeDataViewType();
154  refNodeData_.pyramid_13 = referenceNodeDataViewType();
155  refNodeData_.pyramid_14 = referenceNodeDataViewType();
156 
157  refNodeData_.wedge = referenceNodeDataViewType();
158  refNodeData_.wedge_15 = referenceNodeDataViewType();
159  refNodeData_.wedge_18 = referenceNodeDataViewType();
160  } );
161 
162  isReferenceNodeDataSet_ = true;
163  }
164 
165  template<typename SpT>
166  template<typename cellCenterValueType, class ...cellCenterProperties,
167  typename cellVertexValueType, class ...cellVertexProperties>
168  void
170  getReferenceCellCenter( Kokkos::DynRankView<cellCenterValueType,cellCenterProperties...> cellCenter,
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." );
176 
177  INTREPID2_TEST_FOR_EXCEPTION( cellCenter.rank() != 1, std::invalid_argument,
178  ">>> ERROR (Intrepid2::CellTools::getReferenceCellCenter): cellCenter must have rank 1." );
179 
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." );
182 
183  INTREPID2_TEST_FOR_EXCEPTION( cellVertex.rank() != 1, std::invalid_argument,
184  ">>> ERROR (Intrepid2::CellTools::getReferenceCellCenter): cellVertex must have rank 1." );
185 
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." );
188 #endif
189  const ordinal_type vertexCount = cell.getVertexCount();
190  const ordinal_type dim = cell.getDimension();
191 
192  for (ordinal_type i=0;i<dim;++i) {
193  cellCenter(i) = 0;
194  for (ordinal_type vertOrd=0;vertOrd<vertexCount;++vertOrd) {
195  getReferenceVertex(cellVertex, cell, vertOrd);
196  cellCenter(i) += cellVertex(i);
197  }
198  cellCenter(i) /= vertexCount;
199  }
200  }
201 
202 
203  template<typename SpT>
204  template<typename cellVertexValueType, class ...cellVertexProperties>
205  void
207  getReferenceVertex( Kokkos::DynRankView<cellVertexValueType,cellVertexProperties...> cellVertex,
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." );
213 
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." );
216 
217  INTREPID2_TEST_FOR_EXCEPTION( cellVertex.rank() != 1, std::invalid_argument,
218  ">>> ERROR (Intrepid2::CellTools::getReferenceNode): cellNodes must have rank 1." );
219 
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." );
222 #endif
223  getReferenceNode(cellVertex,
224  cell,
225  vertexOrd);
226  }
227 
228 
229  template<typename SpT>
230  template<typename subcellVertexValueType, class ...subcellVertexProperties>
231  void
233  getReferenceSubcellVertices( Kokkos::DynRankView<subcellVertexValueType,subcellVertexProperties...> subcellVertices,
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." );
240 
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." );
243 
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." );
246 
247  // Verify subcellVertices rank and dimensions
248  INTREPID2_TEST_FOR_EXCEPTION( subcellVertices.rank() != 2, std::invalid_argument,
249  ">>> ERROR (Intrepid2::CellTools::getReferenceSubcellVertices): subcellVertieces must have rank 2." );
250 
251  // need to match to node count as it invokes getReferenceSubcellNodes
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." );
254 
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." );
257 #endif
258  getReferenceSubcellNodes(subcellVertices,
259  subcellDim,
260  subcellOrd,
261  parentCell);
262  }
263 
264 
265  template<typename SpT>
266  template<typename cellNodeValueType, class ...cellNodeProperties>
267  void
269  getReferenceNode( Kokkos::DynRankView<cellNodeValueType,cellNodeProperties...> cellNode,
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." );
275 
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." );
278 
279  INTREPID2_TEST_FOR_EXCEPTION( cellNode.rank() != 1, std::invalid_argument,
280  ">>> ERROR (Intrepid2::CellTools::getReferenceNode): cellNode must have rank 1." );
281 
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." );
284 #endif
285 
286  if (!isReferenceNodeDataSet_)
287  setReferenceNodeData();
288 
289  referenceNodeDataViewType ref;
290 
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;
298 
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;
304 
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;
311 
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;
316 
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;
320 
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;
324 
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;
328 
329  default: {
330  INTREPID2_TEST_FOR_EXCEPTION( true, std::invalid_argument,
331  ">>> ERROR (Intrepid2::CellTools::getReferenceNode): invalid cell topology.");
332  }
333  }
334 
335  // subview version; this is dangerous that users get control over the static data
336  // cellNode = Kokkos::subdynrankview( ref, nodeOrd, Kokkos::ALL() );
337 
338  // let's copy;
339  const ordinal_type dim = cell.getDimension();
340 
341  for (ordinal_type i=0;i<dim;++i)
342  cellNode(i) = ref(nodeOrd, i);
343  }
344 
345  template<typename SpT>
346  template<typename subcellNodeValueType, class ...subcellNodeProperties>
347  void
349  getReferenceSubcellNodes( Kokkos::DynRankView<subcellNodeValueType,subcellNodeProperties...> subcellNodes,
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.");
356 
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.");
359 
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.");
362 
363  // Verify subcellNodes rank and dimensions
364  INTREPID2_TEST_FOR_EXCEPTION( subcellNodes.rank() != 2, std::invalid_argument,
365  ">>> ERROR (Intrepid2::CellTools::getReferenceSubcellNodes): subcellNodes must have rank 2.");
366 
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.");
369 
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.");
372 #endif
373 
374  // Find how many cellWorkset does the specified subcell have.
375  const auto subcNodeCount = parentCell.getNodeCount(subcellDim, subcellOrd);
376 
377  // Loop over subcell cellWorkset
378  for (size_type subcNodeOrd=0;subcNodeOrd<subcNodeCount;++subcNodeOrd) {
379  // Get the node number relative to the parent reference cell
380  const auto cellNodeOrd = parentCell.getNodeMap(subcellDim, subcellOrd, subcNodeOrd);
381 
382  auto dst = Kokkos::subdynrankview(subcellNodes, subcNodeOrd, Kokkos::ALL());
383  getReferenceNode(dst, parentCell, cellNodeOrd);
384  }
385  }
386 
387  template<typename SpT>
388  template<typename refEdgeTangentValueType, class ...refEdgeTangentProperties>
389  void
391  getReferenceEdgeTangent( Kokkos::DynRankView<refEdgeTangentValueType,refEdgeTangentProperties...> refEdgeTangent,
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::getReferenceEdgeTangent): two or three-dimensional parent cell required");
398 
399  INTREPID2_TEST_FOR_EXCEPTION( refEdgeTangent.rank() != 1, std::invalid_argument,
400  ">>> ERROR (Intrepid2::CellTools::getReferenceEdgeTangent): rank = 1 required for output arrays");
401 
402  INTREPID2_TEST_FOR_EXCEPTION( refEdgeTangent.extent(0) != parentCell.getDimension(), std::invalid_argument,
403  ">>> ERROR (Intrepid2::CellTools::getReferenceEdgeTangent): output array size is required to match space dimension");
404 
405  INTREPID2_TEST_FOR_EXCEPTION( edgeOrd < 0 ||
406  edgeOrd >= static_cast<ordinal_type>(parentCell.getSubcellCount(1)), std::invalid_argument,
407  ">>> ERROR (Intrepid2::CellTools::getReferenceEdgeTangent): edge ordinal out of bounds");
408 
409 #endif
410  // Edge parametrizations are computed in setSubcellParametrization and stored in rank-3 array
411  // (subcOrd, coordinate, coefficient)
412  subcellParamViewType edgeMap;
413  getSubcellParametrization(edgeMap, 1, parentCell);
414 
415  // All ref. edge maps have affine coordinate functions: f_dim(u) = C_0(dim) + C_1(dim)*u,
416  // => edge Tangent: -> C_1(*)
417  const ordinal_type dim = parentCell.getDimension();
418  for (ordinal_type i=0;i<dim;++i)
419  refEdgeTangent(i) = edgeMap(edgeOrd, i, 1);
420  }
421 
422 
423  template<typename SpT>
424  template<typename refFaceTanUValueType, class ...refFaceTanUProperties,
425  typename refFaceTanVValueType, class ...refFaceTanVProperties>
426  void
428  getReferenceFaceTangents( Kokkos::DynRankView<refFaceTanUValueType,refFaceTanUProperties...> refFaceTanU,
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");
435 
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");
438 
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");
441 
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");
444 
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");
447 #endif
448 
449  // Face parametrizations are computed in setSubcellParametrization and stored in rank-3 array
450  // (subcOrd, coordinate, coefficient): retrieve this array
451  subcellParamViewType faceMap;
452  getSubcellParametrization(faceMap, 2, parentCell);
453 
454  /* All ref. face maps have affine coordinate functions: f_dim(u,v) = C_0(dim) + C_1(dim)*u + C_2(dim)*v
455  * ` => Tangent vectors are: refFaceTanU -> C_1(*); refFaceTanV -> C_2(*)
456  */
457 
458  // set refFaceTanU -> C_1(*)
459  // set refFaceTanV -> C_2(*)
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);
464  }
465  }
466 
467  template<typename SpT>
468  template<typename refSideNormalValueType, class ...refSideNormalProperties>
469  void
471  getReferenceSideNormal( Kokkos::DynRankView<refSideNormalValueType,refSideNormalProperties...> refSideNormal,
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");
478 
479  // Check side ordinal: by definition side is subcell whose dimension = parentCell.getDimension()-1
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");
482 #endif
483 
484  const auto dim = parentCell.getDimension();
485  if (dim == 2) {
486  // 2D parent cells: side = 1D subcell (edge), call the edge tangent method and rotate tangents
487  getReferenceEdgeTangent(refSideNormal, sideOrd, parentCell);
488 
489  // rotate t(t1, t2) to get n(t2, -t1) so that (n,t) is positively oriented: det(n1,n2/t1,t2)>0
490  refSideNormalValueType tmp = refSideNormal(0);
491  refSideNormal(0) = refSideNormal(1);
492  refSideNormal(1) = -tmp;
493  } else {
494  // 3D parent cell: side = 2D subcell (face), call the face normal method.
495  getReferenceFaceNormal(refSideNormal, sideOrd, parentCell);
496  }
497  }
498 
499 
500  template<typename SpT>
501  template<typename refFaceNormalValueType, class ...refFaceNormalProperties>
502  void
504  getReferenceFaceNormal( Kokkos::DynRankView<refFaceNormalValueType,refFaceNormalProperties...> refFaceNormal,
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");
510 
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");
513 
514  INTREPID2_TEST_FOR_EXCEPTION( refFaceNormal.rank() != 1, std::invalid_argument,
515  ">>> ERROR (Intrepid2::CellTools::getReferenceFaceNormal): rank = 1 required for output array");
516 
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");
519 #endif
520 
521  // Reference face normal = vector product of reference face tangents. Allocate temp FC storage:
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);
528 
529  RealSpaceTools<SpT>::vecprod(refFaceNormal, refFaceTanU, refFaceTanV);
530  }
531 
532  template<typename SpT>
533  template<typename edgeTangentValueType, class ...edgeTangentProperties,
534  typename worksetJacobianValueType, class ...worksetJacobianProperties>
535  void
537  getPhysicalEdgeTangents( Kokkos::DynRankView<edgeTangentValueType,edgeTangentProperties...> edgeTangents,
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." );
545 
546  // (1) edgeTangents is rank-3 (C,P,D) and D=2, or 3 is 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." );
552 
553  // (2) worksetJacobians in rank-4 (C,P,D,D) and D=2, or 3 is required
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." );
561 
562  // (4) cross-check array dimensions: edgeTangents (C,P,D) vs. worksetJacobians (C,P,D,D)
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)." );
566  }
567 #endif
568 
569  // Storage for constant reference edge tangent: rank-1 (D) arrays
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);
575 
576  RealSpaceTools<SpT>::matvec(edgeTangents, worksetJacobians, refEdgeTan);
577  }
578 
579 
580  template<typename SpT>
581  template<typename faceTanUValueType, class ...faceTanUProperties,
582  typename faceTanVValueType, class ...faceTanVProperties,
583  typename worksetJacobianValueType, class ...worksetJacobianProperties>
584  void
586  getPhysicalFaceTangents( Kokkos::DynRankView<faceTanUValueType,faceTanUProperties...> faceTanU,
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");
594 
595  // (1) faceTanU and faceTanV are rank-3 (C,P,D) and D=3 is 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." );
599 
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." );
603 
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." );
607  }
608 
609  // (3) worksetJacobians in rank-4 (C,P,D,D) and D=3 is required
610  INTREPID2_TEST_FOR_EXCEPTION( worksetJacobians.rank() != 4, std::invalid_argument,
611  ">>> ERROR (Intrepid2::CellTools::getPhysicalFaceTangents): worksetJacobians must have rank 4." );
612 
613  INTREPID2_TEST_FOR_EXCEPTION( worksetJacobians.extent(2) != 3, std::invalid_argument,
614  ">>> ERROR (Intrepid2::CellTools::getPhysicalFaceTangents): worksetJacobians dimension(2) must be 3." );
615 
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." );
618 
619  // (4) cross-check array dimensions: faceTanU (C,P,D) vs. worksetJacobians (C,P,D,D)
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." );
623  }
624 #endif
625 
626  // Temp storage for the pair of constant ref. face tangents: rank-1 (D) arrays
627  const auto dim = parentCell.getDimension();
628 
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);
632 
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);
636 
637  getReferenceFaceTangents(refFaceTanU, refFaceTanV, worksetFaceOrd, parentCell);
638 
639  RealSpaceTools<SpT>::matvec(faceTanU, worksetJacobians, refFaceTanU);
640  RealSpaceTools<SpT>::matvec(faceTanV, worksetJacobians, refFaceTanV);
641  }
642 
643 
644  template<typename SpT>
645  template<typename sideNormalValueType, class ...sideNormalProperties,
646  typename worksetJacobianValueType, class ...worksetJacobianProperties>
647  void
649  getPhysicalSideNormals( Kokkos::DynRankView<sideNormalValueType,sideNormalProperties...> sideNormals,
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");
657 
658  // Check side ordinal: by definition side is subcell whose dimension = parentCell.getDimension()-1
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");
662 #endif
663  const auto dim = parentCell.getDimension();
664 
665  if (dim == 2) {
666  // compute edge tangents and rotate it
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);
674 
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;
678 
679  RealSpaceTools<SpT>::matvec(sideNormals, rotation, edgeTangents);
680  } else {
681  getPhysicalFaceNormals(sideNormals, worksetJacobians, worksetSideOrd, parentCell);
682  }
683  }
684 
685 
686  template<typename SpT>
687  template<typename faceNormalValueType, class ...faceNormalProperties,
688  typename worksetJacobianValueType, class ...worksetJacobianProperties>
689  void
691  getPhysicalFaceNormals( Kokkos::DynRankView<faceNormalValueType,faceNormalProperties...> faceNormals,
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." );
698 
699  // (1) faceNormals is rank-3 (C,P,D) and D=3 is 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." );
704 
705  // (3) worksetJacobians in rank-4 (C,P,D,D) and D=3 is required
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)." );
712 
713  // (4) cross-check array dimensions: faceNormals (C,P,D) vs. worksetJacobians (C,P,D,D)
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)." );
717  }
718 #endif
719 
720  // this should be provided from users
721  // Storage for physical face tangents: rank-3 (C,P,D) arrays
722  const auto worksetSize = worksetJacobians.extent(0);
723  const auto facePtCount = worksetJacobians.extent(1);
724  const auto dim = parentCell.getDimension();
725 
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);
730 
731  getPhysicalFaceTangents(faceTanU, faceTanV,
732  worksetJacobians,
733  worksetFaceOrd,
734  parentCell);
735 
736  RealSpaceTools<SpT>::vecprod(faceNormals, faceTanU, faceTanV);
737  }
738 
739 
740  template<typename SpT>
741  bool
743  isReferenceNodeDataSet_ = false;
744 
745  template<typename SpT>
748  refNodeData_ = typename CellTools<SpT>::ReferenceNodeData();
749 
750  template<typename SpT>
753  refNodeDataStatic_ = {
754  // line
755  { // 2
756  {-1.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}
757  },
758  { // 3
759  {-1.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 0.0, 0.0}
760  },
761  // triangle
762  { // 3
763  { 0.0, 0.0, 0.0}, { 1.0, 0.0, 0.0}, { 0.0, 1.0, 0.0}
764  },
765  { // 4
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}
767  },
768  { // 6
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}
771  },
772  // quad
773  { // 4
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}
775  },
776  { // 8
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}
779  },
780  { // 9
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}
783  },
784  // tet
785  { // 4
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}
787  },
788  { // 8
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}
791  },
792  { // 10
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}
795  },
796  { // 11
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}
799  },
800  // hex
801  { // 8
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}
804  },
805  { // 20
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}
811  },
812  { // 27
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},
818  { 0.0, 0.0, 0.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}
820  },
821  // pyramid
822  { // 5
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}
824  },
825  { // 13
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}
829  },
830  { // 14
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}
834  },
835  // wedge
836  { // 6
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}
838  },
839  { // 15
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}
843  },
844  { // 18
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}
849  }
850  };
851 
852 }
853 
854 #endif
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 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 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 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.
Reference node data for each supported topology.
A stateless class for operations on cell data. Provides methods for:
static void matvec(Kokkos::DynRankView< matVecValueType, matVecProperties...> matVecs, const Kokkos::DynRankView< inMatValueType, inMatProperties...> inMats, const Kokkos::DynRankView< inVecValueType, inVecProperties...> inVecs)
Matrix-vector left multiply using multidimensional arrays: matVec = inMat * inVec.
static void getPhysicalFaceTangents(Kokkos::DynRankView< faceTanUValueType, faceTanUProperties...> faceTanU, Kokkos::DynRankView< faceTanVValueType, faceTanVProperties...> faceTanV, const Kokkos::DynRankView< worksetJacobianValueType, worksetJacobianProperties...> worksetJacobians, const ordinal_type worksetFaceOrd, const shards::CellTopology parentCell)
Computes non-normalized tangent vector pairs to physical faces in a face workset ; (see Subcell works...
static void 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 vecprod(Kokkos::DynRankView< vecProdValueType, vecProdProperties...> vecProd, const Kokkos::DynRankView< inLeftValueType, inLeftProperties...> inLeft, const Kokkos::DynRankView< inRightValueType, inRightProperties...> inRight)
Vector product using multidimensional arrays: vecProd = inVecLeft x inVecRight
static void setReferenceNodeData()
Set reference node coordinates for supported topologies.
static void getReferenceFaceTangents(Kokkos::DynRankView< refFaceTanUValueType, refFaceTanUProperties...> refFaceTanU, Kokkos::DynRankView< refFaceTanVValueType, refFaceTanVProperties...> refFaceTanV, const ordinal_type faceOrd, const shards::CellTopology parentCell)
Computes pairs of constant tangent vectors to faces of a 3D reference cells.
static void getPhysicalSideNormals(Kokkos::DynRankView< sideNormalValueType, sideNormalProperties...> sideNormals, const Kokkos::DynRankView< worksetJacobianValueType, worksetJacobianProperties...> worksetJacobians, const ordinal_type worksetSideOrd, const shards::CellTopology parentCell)
Computes non-normalized normal vectors to physical sides in a side workset .
static void getReferenceCellCenter(Kokkos::DynRankView< cellCenterValueType, cellCenterProperties...> cellCenter, Kokkos::DynRankView< cellVertexValueType, cellVertexProperties...> cellVertex, const shards::CellTopology cell)
Computes the Cartesian coordinates of reference cell center.
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...
Reference node containers for each supported topology.
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 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 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.