Intrepid2
Intrepid2_PointTools.hpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Intrepid2 Package
4 //
5 // Copyright 2007 NTESS and the Intrepid2 contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
17 #ifndef __INTREPID2_POINTTOOLS_HPP__
18 #define __INTREPID2_POINTTOOLS_HPP__
19 
20 #include "Intrepid2_ConfigDefs.hpp"
21 
22 #include "Intrepid2_Types.hpp"
23 #include "Intrepid2_Utils.hpp"
24 
25 #include "Shards_CellTopology.hpp"
26 
27 #include "Intrepid2_Polylib.hpp"
28 #include "Intrepid2_CellTools.hpp"
29 
30 #include "Kokkos_Core.hpp"
31 
32 namespace Intrepid2 {
33 
166 class PointTools {
167 public:
185  inline
186  static ordinal_type
187  getLatticeSize( const shards::CellTopology cellType,
188  const ordinal_type order,
189  const ordinal_type offset = 0 );
190 
207  template<typename pointValueType, class ...pointProperties>
208  static void
209  getLattice( Kokkos::DynRankView<pointValueType,pointProperties...> points,
210  const shards::CellTopology cellType,
211  const ordinal_type order,
212  const ordinal_type offset = 0 ,
213  const EPointType pointType = POINTTYPE_EQUISPACED );
214 
229  template<typename pointValueType, class ...pointProperties>
230  static void
231  getLatticeLine( Kokkos::DynRankView<pointValueType,pointProperties...> points,
232  const ordinal_type order,
233  const ordinal_type offset = 0 ,
234  const EPointType pointType = POINTTYPE_EQUISPACED );
235 
250  template<typename pointValueType, class ...pointProperties>
251  static void
252  getLatticeTriangle( Kokkos::DynRankView<pointValueType,pointProperties...> points,
253  const ordinal_type order,
254  const ordinal_type offset = 0 ,
255  const EPointType pointType = POINTTYPE_EQUISPACED );
256 
271  template<typename pointValueType, class ...pointProperties>
272  static void
273  getLatticeTetrahedron( Kokkos::DynRankView<pointValueType,pointProperties...> points,
274  const ordinal_type order,
275  const ordinal_type offset = 0 ,
276  const EPointType pointType = POINTTYPE_EQUISPACED );
277 
290  template<typename pointValueType, class ...pointProperties>
291  static void
292  getLatticePyramid( Kokkos::DynRankView<pointValueType,pointProperties...> points,
293  const ordinal_type order,
294  const ordinal_type offset = 0 ,
295  const EPointType pointType = POINTTYPE_EQUISPACED );
296 
302  template<typename pointValueType, class ...pointProperties>
303  static void getGaussPoints( Kokkos::DynRankView<pointValueType,pointProperties...> points,
304  const ordinal_type order );
305 
306 private:
307 
322  template<typename pointValueType, class ...pointProperties>
323  static void
324  getEquispacedLatticeLine( Kokkos::DynRankView<pointValueType,pointProperties...> points,
325  const ordinal_type order,
326  const ordinal_type offset = 0 );
327 
328 
343  template<typename pointValueType, class ...pointProperties>
344  static void
345  getWarpBlendLatticeLine( Kokkos::DynRankView<pointValueType,pointProperties...> points,
346  const ordinal_type order ,
347  const ordinal_type offset = 0 );
348 
349 
350  // /* \brief Converts Cartesian coordinates to barycentric coordinates
351  // on a batch of triangles.
352  // The input array cartValues is (C,P,2)
353  // The output array baryValues is (C,P,3).
354  // The input array vertices is (C,3,2), where
355  // \code
356  // C - num. integration domains
357  // P - number of points per cell
358  // \endcode
359 
360  // \param baryValues [out] - Output array of barycentric coords
361  // \param cartValues [in] - Input array of Cartesian coords
362  // \param vertices [out] - Vertices of each cell.
363 
364  // */
365  // template<class Scalar, class ArrayTypeOut, class ArrayTypeIn1, class ArrayTypeIn2>
366  // static void cartToBaryTriangle( ArrayTypeOut & baryValues ,
367  // const ArrayTypeIn1 & cartValues ,
368  // const ArrayTypeIn2 & vertices );
369 
370 
371  // /* \brief Converts barycentric coordinates to Cartesian coordinates
372  // on a batch of triangles.
373  // The input array baryValues is (C,P,3)
374  // The output array cartValues is (C,P,2).
375  // The input array vertices is (C,3,2), where
376  // \code
377  // C - num. integration domains
378  // P - number of points per cell
379  // D - is the spatial dimension
380  // \endcode
381 
382  // \param baryValues [out] - Output array of barycentric coords
383  // \param cartValues [in] - Input array of Cartesian coords
384  // \param vertices [out] - Vertices of each cell.
385 
386  // */
387  // template<class Scalar, class ArrayTypeOut, class ArrayTypeIn1, class ArrayTypeIn2>
388  // static void baryToCartTriangle( ArrayTypeOut & cartValues ,
389  // const ArrayTypeIn1 & baryValues ,
390  // const ArrayTypeIn2 & vertices );
391 
392 
393  // /* \brief Converts Cartesian coordinates to barycentric coordinates
394  // on a batch of tetrahedra.
395  // The input array cartValues is (C,P,3)
396  // The output array baryValues is (C,P,4).
397  // The input array vertices is (C,4,3), where
398  // \code
399  // C - num. integration domains
400  // P - number of points per cell
401  // D - is the spatial dimension
402  // \endcode
403 
404  // \param baryValues [out] - Output array of barycentric coords
405  // \param cartValues [in] - Input array of Cartesian coords
406  // \param vertices [out] - Vertices of each cell.
407 
408  // */
409  // template<class Scalar, class ArrayTypeOut, class ArrayTypeIn1, class ArrayTypeIn2>
410  // static void cartToBaryTetrahedron( ArrayTypeOut & baryValues ,
411  // const ArrayTypeIn1 & cartValues ,
412  // const ArrayTypeIn2 & vertices );
413 
414  // /* \brief Converts barycentric coordinates to Cartesian coordinates
415  // on a batch of tetrahedra.
416  // The input array baryValues is (C,P,4)
417  // The output array cartValues is (C,P,3).
418  // The input array vertices is (C,4,3), where
419  // \code
420  // C - num. integration domains
421  // P - number of points per cell
422  // D - is the spatial dimension
423  // \endcode
424 
425  // \param baryValues [out] - Output array of barycentric coords
426  // \param cartValues [in] - Input array of Cartesian coords
427  // \param vertices [out] - Vertices of each cell.
428 
429  // */
430  // template<class Scalar, class ArrayTypeOut, class ArrayTypeIn1, class ArrayTypeIn2>
431  // static void baryToCartTetrahedron( ArrayTypeOut & cartValues ,
432  // const ArrayTypeIn1 & baryValues ,
433  // const ArrayTypeIn2 & vertices );
434 
435 
436 
437 
452  template<typename pointValueType, class ...pointProperties>
453  static void
454  getEquispacedLatticeTriangle( Kokkos::DynRankView<pointValueType,pointProperties...> points,
455  const ordinal_type order,
456  const ordinal_type offset = 0 );
457 
470  template<typename pointValueType, class ...pointProperties>
471  static void
472  getEquispacedLatticeTetrahedron( Kokkos::DynRankView<pointValueType,pointProperties...> points,
473  const ordinal_type order ,
474  const ordinal_type offset = 0 );
475 
476 
489  template<typename pointValueType, class ...pointProperties>
490  static void
491  getEquispacedLatticePyramid( Kokkos::DynRankView<pointValueType,pointProperties...> points,
492  const ordinal_type order,
493  const ordinal_type offset = 0 );
494 
495 
502  template<typename pointValueType, class ...pointProperties>
503  static void
504  warpFactor( Kokkos::DynRankView<pointValueType,pointProperties...> warp ,
505  const ordinal_type order ,
506  const Kokkos::DynRankView<pointValueType,pointProperties...> xnodes ,
507  const Kokkos::DynRankView<pointValueType,pointProperties...> xout
508  );
509 
524  template<typename pointValueType, class ...pointProperties>
525  static void
526  getWarpBlendLatticeTriangle( Kokkos::DynRankView<pointValueType,pointProperties...> points,
527  const ordinal_type order ,
528  const ordinal_type offset = 0 );
529 
543  template<typename pointValueType, class ...pointProperties>
544  static void
545  getWarpBlendLatticeTetrahedron( Kokkos::DynRankView<pointValueType,pointProperties...> points,
546  const ordinal_type order ,
547  const ordinal_type offset = 0 );
548 
549 
560  template<typename pointValueType, class ...pointProperties>
561  static void
562  warpShiftFace3D( Kokkos::DynRankView<pointValueType,pointProperties...> dxy,
563  const ordinal_type order ,
564  const pointValueType pval ,
565  const Kokkos::DynRankView<pointValueType,pointProperties...> L1,
566  const Kokkos::DynRankView<pointValueType,pointProperties...> L2,
567  const Kokkos::DynRankView<pointValueType,pointProperties...> L3,
568  const Kokkos::DynRankView<pointValueType,pointProperties...> L4
569  );
570 
580  template<typename pointValueType, class ...pointProperties>
581  static void
582  evalshift( Kokkos::DynRankView<pointValueType,pointProperties...> dxy ,
583  const ordinal_type order ,
584  const pointValueType pval ,
585  const Kokkos::DynRankView<pointValueType,pointProperties...> L1 ,
586  const Kokkos::DynRankView<pointValueType,pointProperties...> L2 ,
587  const Kokkos::DynRankView<pointValueType,pointProperties...> L3
588  );
589 
596  template<typename pointValueType, class ...pointProperties>
597  static void
598  evalwarp( Kokkos::DynRankView<pointValueType,pointProperties...> warp ,
599  const ordinal_type order ,
600  const Kokkos::DynRankView<pointValueType,pointProperties...> xnodes ,
601  const Kokkos::DynRankView<pointValueType,pointProperties...> xout );
602 };
603 
604 }
605 
607 
608 #endif
static void evalwarp(Kokkos::DynRankView< pointValueType, pointProperties...> warp, const ordinal_type order, const Kokkos::DynRankView< pointValueType, pointProperties...> xnodes, const Kokkos::DynRankView< pointValueType, pointProperties...> xout)
Used internally to compute the warp on edges of a triangle in warp-blend points.
static void evalshift(Kokkos::DynRankView< pointValueType, pointProperties...> dxy, const ordinal_type order, const pointValueType pval, const Kokkos::DynRankView< pointValueType, pointProperties...> L1, const Kokkos::DynRankView< pointValueType, pointProperties...> L2, const Kokkos::DynRankView< pointValueType, pointProperties...> L3)
Used internally to evaluate the point shift for warp-blend points on faces of tets.
Header file for Intrepid2::Polylib class providing orthogonal polynomial calculus and interpolation...
static void warpFactor(Kokkos::DynRankView< pointValueType, pointProperties...> warp, const ordinal_type order, const Kokkos::DynRankView< pointValueType, pointProperties...> xnodes, const Kokkos::DynRankView< pointValueType, pointProperties...> xout)
interpolates Warburton&#39;s warp function on the line
Header function for Intrepid2::Util class and other utility functions.
static void getLatticeTetrahedron(Kokkos::DynRankView< pointValueType, pointProperties...> points, const ordinal_type order, const ordinal_type offset=0, const EPointType pointType=POINTTYPE_EQUISPACED)
Computes a lattice of points of a given order on a reference tetrahedron. The output array is (P...
static void getLattice(Kokkos::DynRankView< pointValueType, pointProperties...> points, const shards::CellTopology cellType, const ordinal_type order, const ordinal_type offset=0, const EPointType pointType=POINTTYPE_EQUISPACED)
Computes a lattice of points of a given order on a reference simplex, quadrilateral or hexahedron (cu...
static void getEquispacedLatticeTetrahedron(Kokkos::DynRankView< pointValueType, pointProperties...> points, const ordinal_type order, const ordinal_type offset=0)
Computes an equispaced lattice of a given order on the reference tetrahedron. The output array is (P...
static void getLatticePyramid(Kokkos::DynRankView< pointValueType, pointProperties...> points, const ordinal_type order, const ordinal_type offset=0, const EPointType pointType=POINTTYPE_EQUISPACED)
Computes a lattice of points of a given order on a reference pyramid. The output array is (P...
static void getWarpBlendLatticeLine(Kokkos::DynRankView< pointValueType, pointProperties...> points, const ordinal_type order, const ordinal_type offset=0)
Returns the Gauss-Lobatto points of a given order on the reference line [-1,1]. The output array is (...
Contains definitions of custom data types in Intrepid2.
static void getWarpBlendLatticeTriangle(Kokkos::DynRankView< pointValueType, pointProperties...> points, const ordinal_type order, const ordinal_type offset=0)
Returns Warburton&#39;s warp-blend points of a given order on the reference triangle. The output array is...
static void getEquispacedLatticeLine(Kokkos::DynRankView< pointValueType, pointProperties...> points, const ordinal_type order, const ordinal_type offset=0)
Computes an equispaced lattice of a given order on the reference line [-1,1]. The output array is (P...
static void getEquispacedLatticePyramid(Kokkos::DynRankView< pointValueType, pointProperties...> points, const ordinal_type order, const ordinal_type offset=0)
Computes an equispaced lattice of a given order on the reference pyramid. The output array is (P...
static void getLatticeTriangle(Kokkos::DynRankView< pointValueType, pointProperties...> points, const ordinal_type order, const ordinal_type offset=0, const EPointType pointType=POINTTYPE_EQUISPACED)
Computes a lattice of points of a given order on a reference triangle. The output array is (P...
Utility class that provides methods for calculating distributions of points on different cells...
static void getGaussPoints(Kokkos::DynRankView< pointValueType, pointProperties...> points, const ordinal_type order)
static void warpShiftFace3D(Kokkos::DynRankView< pointValueType, pointProperties...> dxy, const ordinal_type order, const pointValueType pval, const Kokkos::DynRankView< pointValueType, pointProperties...> L1, const Kokkos::DynRankView< pointValueType, pointProperties...> L2, const Kokkos::DynRankView< pointValueType, pointProperties...> L3, const Kokkos::DynRankView< pointValueType, pointProperties...> L4)
This is used internally to compute the tetrahedral warp-blend points one each face.
static void getLatticeLine(Kokkos::DynRankView< pointValueType, pointProperties...> points, const ordinal_type order, const ordinal_type offset=0, const EPointType pointType=POINTTYPE_EQUISPACED)
Computes a lattice of points of a given order on a reference line. The output array is (P...
static void getWarpBlendLatticeTetrahedron(Kokkos::DynRankView< pointValueType, pointProperties...> points, const ordinal_type order, const ordinal_type offset=0)
Returns Warburton&#39;s warp-blend points of a given order on the reference tetrahedron. The output array is (P,3), where.
Header file for the Intrepid2::CellTools class.
static ordinal_type getLatticeSize(const shards::CellTopology cellType, const ordinal_type order, const ordinal_type offset=0)
Computes the number of points in a lattice of a given order on a simplex (currently disabled for othe...
Definition file for point tool utilities for barycentric coordinates and lattices.
static void getEquispacedLatticeTriangle(Kokkos::DynRankView< pointValueType, pointProperties...> points, const ordinal_type order, const ordinal_type offset=0)
Computes an equispaced lattice of a given order on the reference triangle. The output array is (P...