Intrepid2
Intrepid2_PointTools.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 
50 #ifndef __INTREPID2_POINTTOOLS_HPP__
51 #define __INTREPID2_POINTTOOLS_HPP__
52 
53 #include "Intrepid2_ConfigDefs.hpp"
54 
55 #include "Intrepid2_Types.hpp"
56 #include "Intrepid2_Utils.hpp"
57 
58 #include "Shards_CellTopology.hpp"
59 
60 #include "Intrepid2_Polylib.hpp"
61 #include "Intrepid2_CellTools.hpp"
62 
63 #include "Kokkos_Core.hpp"
64 
65 namespace Intrepid2 {
66 
199 class PointTools {
200 public:
218  inline
219  static ordinal_type
220  getLatticeSize( const shards::CellTopology cellType,
221  const ordinal_type order,
222  const ordinal_type offset = 0 );
223 
240  template<typename pointValueType, class ...pointProperties>
241  static void
242  getLattice( Kokkos::DynRankView<pointValueType,pointProperties...> points,
243  const shards::CellTopology cellType,
244  const ordinal_type order,
245  const ordinal_type offset = 0 ,
246  const EPointType pointType = POINTTYPE_EQUISPACED );
247 
262  template<typename pointValueType, class ...pointProperties>
263  static void
264  getLatticeLine( Kokkos::DynRankView<pointValueType,pointProperties...> points,
265  const ordinal_type order,
266  const ordinal_type offset = 0 ,
267  const EPointType pointType = POINTTYPE_EQUISPACED );
268 
283  template<typename pointValueType, class ...pointProperties>
284  static void
285  getLatticeTriangle( Kokkos::DynRankView<pointValueType,pointProperties...> points,
286  const ordinal_type order,
287  const ordinal_type offset = 0 ,
288  const EPointType pointType = POINTTYPE_EQUISPACED );
289 
304  template<typename pointValueType, class ...pointProperties>
305  static void
306  getLatticeTetrahedron( Kokkos::DynRankView<pointValueType,pointProperties...> points,
307  const ordinal_type order,
308  const ordinal_type offset = 0 ,
309  const EPointType pointType = POINTTYPE_EQUISPACED );
310 
323  template<typename pointValueType, class ...pointProperties>
324  static void
325  getLatticePyramid( Kokkos::DynRankView<pointValueType,pointProperties...> points,
326  const ordinal_type order,
327  const ordinal_type offset = 0 ,
328  const EPointType pointType = POINTTYPE_EQUISPACED );
329 
335  template<typename pointValueType, class ...pointProperties>
336  static void getGaussPoints( Kokkos::DynRankView<pointValueType,pointProperties...> points,
337  const ordinal_type order );
338 
339 private:
340 
355  template<typename pointValueType, class ...pointProperties>
356  static void
357  getEquispacedLatticeLine( Kokkos::DynRankView<pointValueType,pointProperties...> points,
358  const ordinal_type order,
359  const ordinal_type offset = 0 );
360 
361 
376  template<typename pointValueType, class ...pointProperties>
377  static void
378  getWarpBlendLatticeLine( Kokkos::DynRankView<pointValueType,pointProperties...> points,
379  const ordinal_type order ,
380  const ordinal_type offset = 0 );
381 
382 
383  // /* \brief Converts Cartesian coordinates to barycentric coordinates
384  // on a batch of triangles.
385  // The input array cartValues is (C,P,2)
386  // The output array baryValues is (C,P,3).
387  // The input array vertices is (C,3,2), where
388  // \code
389  // C - num. integration domains
390  // P - number of points per cell
391  // \endcode
392 
393  // \param baryValues [out] - Output array of barycentric coords
394  // \param cartValues [in] - Input array of Cartesian coords
395  // \param vertices [out] - Vertices of each cell.
396 
397  // */
398  // template<class Scalar, class ArrayTypeOut, class ArrayTypeIn1, class ArrayTypeIn2>
399  // static void cartToBaryTriangle( ArrayTypeOut & baryValues ,
400  // const ArrayTypeIn1 & cartValues ,
401  // const ArrayTypeIn2 & vertices );
402 
403 
404  // /* \brief Converts barycentric coordinates to Cartesian coordinates
405  // on a batch of triangles.
406  // The input array baryValues is (C,P,3)
407  // The output array cartValues is (C,P,2).
408  // The input array vertices is (C,3,2), where
409  // \code
410  // C - num. integration domains
411  // P - number of points per cell
412  // D - is the spatial dimension
413  // \endcode
414 
415  // \param baryValues [out] - Output array of barycentric coords
416  // \param cartValues [in] - Input array of Cartesian coords
417  // \param vertices [out] - Vertices of each cell.
418 
419  // */
420  // template<class Scalar, class ArrayTypeOut, class ArrayTypeIn1, class ArrayTypeIn2>
421  // static void baryToCartTriangle( ArrayTypeOut & cartValues ,
422  // const ArrayTypeIn1 & baryValues ,
423  // const ArrayTypeIn2 & vertices );
424 
425 
426  // /* \brief Converts Cartesian coordinates to barycentric coordinates
427  // on a batch of tetrahedra.
428  // The input array cartValues is (C,P,3)
429  // The output array baryValues is (C,P,4).
430  // The input array vertices is (C,4,3), where
431  // \code
432  // C - num. integration domains
433  // P - number of points per cell
434  // D - is the spatial dimension
435  // \endcode
436 
437  // \param baryValues [out] - Output array of barycentric coords
438  // \param cartValues [in] - Input array of Cartesian coords
439  // \param vertices [out] - Vertices of each cell.
440 
441  // */
442  // template<class Scalar, class ArrayTypeOut, class ArrayTypeIn1, class ArrayTypeIn2>
443  // static void cartToBaryTetrahedron( ArrayTypeOut & baryValues ,
444  // const ArrayTypeIn1 & cartValues ,
445  // const ArrayTypeIn2 & vertices );
446 
447  // /* \brief Converts barycentric coordinates to Cartesian coordinates
448  // on a batch of tetrahedra.
449  // The input array baryValues is (C,P,4)
450  // The output array cartValues is (C,P,3).
451  // The input array vertices is (C,4,3), where
452  // \code
453  // C - num. integration domains
454  // P - number of points per cell
455  // D - is the spatial dimension
456  // \endcode
457 
458  // \param baryValues [out] - Output array of barycentric coords
459  // \param cartValues [in] - Input array of Cartesian coords
460  // \param vertices [out] - Vertices of each cell.
461 
462  // */
463  // template<class Scalar, class ArrayTypeOut, class ArrayTypeIn1, class ArrayTypeIn2>
464  // static void baryToCartTetrahedron( ArrayTypeOut & cartValues ,
465  // const ArrayTypeIn1 & baryValues ,
466  // const ArrayTypeIn2 & vertices );
467 
468 
469 
470 
485  template<typename pointValueType, class ...pointProperties>
486  static void
487  getEquispacedLatticeTriangle( Kokkos::DynRankView<pointValueType,pointProperties...> points,
488  const ordinal_type order,
489  const ordinal_type offset = 0 );
490 
503  template<typename pointValueType, class ...pointProperties>
504  static void
505  getEquispacedLatticeTetrahedron( Kokkos::DynRankView<pointValueType,pointProperties...> points,
506  const ordinal_type order ,
507  const ordinal_type offset = 0 );
508 
509 
522  template<typename pointValueType, class ...pointProperties>
523  static void
524  getEquispacedLatticePyramid( Kokkos::DynRankView<pointValueType,pointProperties...> points,
525  const ordinal_type order,
526  const ordinal_type offset = 0 );
527 
528 
535  template<typename pointValueType, class ...pointProperties>
536  static void
537  warpFactor( Kokkos::DynRankView<pointValueType,pointProperties...> warp ,
538  const ordinal_type order ,
539  const Kokkos::DynRankView<pointValueType,pointProperties...> xnodes ,
540  const Kokkos::DynRankView<pointValueType,pointProperties...> xout
541  );
542 
557  template<typename pointValueType, class ...pointProperties>
558  static void
559  getWarpBlendLatticeTriangle( Kokkos::DynRankView<pointValueType,pointProperties...> points,
560  const ordinal_type order ,
561  const ordinal_type offset = 0 );
562 
576  template<typename pointValueType, class ...pointProperties>
577  static void
578  getWarpBlendLatticeTetrahedron( Kokkos::DynRankView<pointValueType,pointProperties...> points,
579  const ordinal_type order ,
580  const ordinal_type offset = 0 );
581 
582 
593  template<typename pointValueType, class ...pointProperties>
594  static void
595  warpShiftFace3D( Kokkos::DynRankView<pointValueType,pointProperties...> dxy,
596  const ordinal_type order ,
597  const pointValueType pval ,
598  const Kokkos::DynRankView<pointValueType,pointProperties...> L1,
599  const Kokkos::DynRankView<pointValueType,pointProperties...> L2,
600  const Kokkos::DynRankView<pointValueType,pointProperties...> L3,
601  const Kokkos::DynRankView<pointValueType,pointProperties...> L4
602  );
603 
613  template<typename pointValueType, class ...pointProperties>
614  static void
615  evalshift( Kokkos::DynRankView<pointValueType,pointProperties...> dxy ,
616  const ordinal_type order ,
617  const pointValueType pval ,
618  const Kokkos::DynRankView<pointValueType,pointProperties...> L1 ,
619  const Kokkos::DynRankView<pointValueType,pointProperties...> L2 ,
620  const Kokkos::DynRankView<pointValueType,pointProperties...> L3
621  );
622 
629  template<typename pointValueType, class ...pointProperties>
630  static void
631  evalwarp( Kokkos::DynRankView<pointValueType,pointProperties...> warp ,
632  const ordinal_type order ,
633  const Kokkos::DynRankView<pointValueType,pointProperties...> xnodes ,
634  const Kokkos::DynRankView<pointValueType,pointProperties...> xout );
635 };
636 
637 }
638 
640 
641 #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...