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 
316  template<typename pointValueType, class ...pointProperties>
317  static void getGaussPoints( Kokkos::DynRankView<pointValueType,pointProperties...> points,
318  const ordinal_type order );
319 
320 private:
321 
336  template<typename pointValueType, class ...pointProperties>
337  static void
338  getEquispacedLatticeLine( Kokkos::DynRankView<pointValueType,pointProperties...> points,
339  const ordinal_type order,
340  const ordinal_type offset = 0 );
341 
342 
357  template<typename pointValueType, class ...pointProperties>
358  static void
359  getWarpBlendLatticeLine( Kokkos::DynRankView<pointValueType,pointProperties...> points,
360  const ordinal_type order ,
361  const ordinal_type offset = 0 );
362 
363 
364  // /* \brief Converts Cartesian coordinates to barycentric coordinates
365  // on a batch of triangles.
366  // The input array cartValues is (C,P,2)
367  // The output array baryValues is (C,P,3).
368  // The input array vertices is (C,3,2), where
369  // \code
370  // C - num. integration domains
371  // P - number of points per cell
372  // \endcode
373 
374  // \param baryValues [out] - Output array of barycentric coords
375  // \param cartValues [in] - Input array of Cartesian coords
376  // \param vertices [out] - Vertices of each cell.
377 
378  // */
379  // template<class Scalar, class ArrayTypeOut, class ArrayTypeIn1, class ArrayTypeIn2>
380  // static void cartToBaryTriangle( ArrayTypeOut & baryValues ,
381  // const ArrayTypeIn1 & cartValues ,
382  // const ArrayTypeIn2 & vertices );
383 
384 
385  // /* \brief Converts barycentric coordinates to Cartesian coordinates
386  // on a batch of triangles.
387  // The input array baryValues is (C,P,3)
388  // The output array cartValues is (C,P,2).
389  // The input array vertices is (C,3,2), where
390  // \code
391  // C - num. integration domains
392  // P - number of points per cell
393  // D - is the spatial dimension
394  // \endcode
395 
396  // \param baryValues [out] - Output array of barycentric coords
397  // \param cartValues [in] - Input array of Cartesian coords
398  // \param vertices [out] - Vertices of each cell.
399 
400  // */
401  // template<class Scalar, class ArrayTypeOut, class ArrayTypeIn1, class ArrayTypeIn2>
402  // static void baryToCartTriangle( ArrayTypeOut & cartValues ,
403  // const ArrayTypeIn1 & baryValues ,
404  // const ArrayTypeIn2 & vertices );
405 
406 
407  // /* \brief Converts Cartesian coordinates to barycentric coordinates
408  // on a batch of tetrahedra.
409  // The input array cartValues is (C,P,3)
410  // The output array baryValues is (C,P,4).
411  // The input array vertices is (C,4,3), where
412  // \code
413  // C - num. integration domains
414  // P - number of points per cell
415  // D - is the spatial dimension
416  // \endcode
417 
418  // \param baryValues [out] - Output array of barycentric coords
419  // \param cartValues [in] - Input array of Cartesian coords
420  // \param vertices [out] - Vertices of each cell.
421 
422  // */
423  // template<class Scalar, class ArrayTypeOut, class ArrayTypeIn1, class ArrayTypeIn2>
424  // static void cartToBaryTetrahedron( ArrayTypeOut & baryValues ,
425  // const ArrayTypeIn1 & cartValues ,
426  // const ArrayTypeIn2 & vertices );
427 
428  // /* \brief Converts barycentric coordinates to Cartesian coordinates
429  // on a batch of tetrahedra.
430  // The input array baryValues is (C,P,4)
431  // The output array cartValues is (C,P,3).
432  // The input array vertices is (C,4,3), where
433  // \code
434  // C - num. integration domains
435  // P - number of points per cell
436  // D - is the spatial dimension
437  // \endcode
438 
439  // \param baryValues [out] - Output array of barycentric coords
440  // \param cartValues [in] - Input array of Cartesian coords
441  // \param vertices [out] - Vertices of each cell.
442 
443  // */
444  // template<class Scalar, class ArrayTypeOut, class ArrayTypeIn1, class ArrayTypeIn2>
445  // static void baryToCartTetrahedron( ArrayTypeOut & cartValues ,
446  // const ArrayTypeIn1 & baryValues ,
447  // const ArrayTypeIn2 & vertices );
448 
449 
450 
451 
466  template<typename pointValueType, class ...pointProperties>
467  static void
468  getEquispacedLatticeTriangle( Kokkos::DynRankView<pointValueType,pointProperties...> points,
469  const ordinal_type order,
470  const ordinal_type offset = 0 );
471 
485  template<typename pointValueType, class ...pointProperties>
486  static void
487  getEquispacedLatticeTetrahedron( Kokkos::DynRankView<pointValueType,pointProperties...> points,
488  const ordinal_type order ,
489  const ordinal_type offset = 0 );
490 
491 
498  template<typename pointValueType, class ...pointProperties>
499  static void
500  warpFactor( Kokkos::DynRankView<pointValueType,pointProperties...> warp ,
501  const ordinal_type order ,
502  const Kokkos::DynRankView<pointValueType,pointProperties...> xnodes ,
503  const Kokkos::DynRankView<pointValueType,pointProperties...> xout
504  );
505 
520  template<typename pointValueType, class ...pointProperties>
521  static void
522  getWarpBlendLatticeTriangle( Kokkos::DynRankView<pointValueType,pointProperties...> points,
523  const ordinal_type order ,
524  const ordinal_type offset = 0 );
525 
539  template<typename pointValueType, class ...pointProperties>
540  static void
541  getWarpBlendLatticeTetrahedron( Kokkos::DynRankView<pointValueType,pointProperties...> points,
542  const ordinal_type order ,
543  const ordinal_type offset = 0 );
544 
545 
556  template<typename pointValueType, class ...pointProperties>
557  static void
558  warpShiftFace3D( Kokkos::DynRankView<pointValueType,pointProperties...> dxy,
559  const ordinal_type order ,
560  const pointValueType pval ,
561  const Kokkos::DynRankView<pointValueType,pointProperties...> L1,
562  const Kokkos::DynRankView<pointValueType,pointProperties...> L2,
563  const Kokkos::DynRankView<pointValueType,pointProperties...> L3,
564  const Kokkos::DynRankView<pointValueType,pointProperties...> L4
565  );
566 
576  template<typename pointValueType, class ...pointProperties>
577  static void
578  evalshift( Kokkos::DynRankView<pointValueType,pointProperties...> dxy ,
579  const ordinal_type order ,
580  const pointValueType pval ,
581  const Kokkos::DynRankView<pointValueType,pointProperties...> L1 ,
582  const Kokkos::DynRankView<pointValueType,pointProperties...> L2 ,
583  const Kokkos::DynRankView<pointValueType,pointProperties...> L3
584  );
585 
592  template<typename pointValueType, class ...pointProperties>
593  static void
594  evalwarp( Kokkos::DynRankView<pointValueType,pointProperties...> warp ,
595  const ordinal_type order ,
596  const Kokkos::DynRankView<pointValueType,pointProperties...> xnodes ,
597  const Kokkos::DynRankView<pointValueType,pointProperties...> xout );
598 };
599 
600 }
601 
603 
604 #endif
Header file for the Intrepid2::CellTools class.
Definition file for point tool utilities for barycentric coordinates and lattices.
Header file for Intrepid2::Polylib class providing orthogonal polynomial calculus and interpolation.
Contains definitions of custom data types in Intrepid2.
Header function for Intrepid2::Util class and other utility functions.
Utility class that provides methods for calculating distributions of points on different cells.
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 (...
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,...
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.
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 warpFactor(Kokkos::DynRankView< pointValueType, pointProperties... > warp, const ordinal_type order, const Kokkos::DynRankView< pointValueType, pointProperties... > xnodes, const Kokkos::DynRankView< pointValueType, pointProperties... > xout)
interpolates Warburton's warp function on the line
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...
static void getWarpBlendLatticeTetrahedron(Kokkos::DynRankView< pointValueType, pointProperties... > points, const ordinal_type order, const ordinal_type offset=0)
Returns Warburton's warp-blend points of a given order on the reference tetrahedron....
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,...
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 getGaussPoints(Kokkos::DynRankView< pointValueType, pointProperties... > points, const ordinal_type order)
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 getWarpBlendLatticeTriangle(Kokkos::DynRankView< pointValueType, pointProperties... > points, const ordinal_type order, const ordinal_type offset=0)
Returns Warburton's warp-blend points of a given order on the reference triangle. The output array is...
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 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 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 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,...