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 
253  template<typename pointValueType, class ...pointProperties>
254  static void getGaussPoints( Kokkos::DynRankView<pointValueType,pointProperties...> points,
255  const ordinal_type order );
256 
257 private:
258 
273  template<typename pointValueType, class ...pointProperties>
274  static void
275  getEquispacedLattice( Kokkos::DynRankView<pointValueType,pointProperties...> points,
276  const shards::CellTopology cellType,
277  const ordinal_type order,
278  const ordinal_type offset = 0 );
279 
280 
295  template<typename pointValueType, class ...pointProperties>
296  static void
297  getWarpBlendLattice( Kokkos::DynRankView<pointValueType,pointProperties...> points,
298  const shards::CellTopology cellType,
299  const ordinal_type order,
300  const ordinal_type offset = 0 );
301 
302 
317  template<typename pointValueType, class ...pointProperties>
318  static void
319  getEquispacedLatticeLine( Kokkos::DynRankView<pointValueType,pointProperties...> points,
320  const ordinal_type order,
321  const ordinal_type offset = 0 );
322 
323 
338  template<typename pointValueType, class ...pointProperties>
339  static void
340  getWarpBlendLatticeLine( Kokkos::DynRankView<pointValueType,pointProperties...> points,
341  const ordinal_type order ,
342  const ordinal_type offset = 0 );
343 
344 
345  // /* \brief Converts Cartesian coordinates to barycentric coordinates
346  // on a batch of triangles.
347  // The input array cartValues is (C,P,2)
348  // The output array baryValues is (C,P,3).
349  // The input array vertices is (C,3,2), where
350  // \code
351  // C - num. integration domains
352  // P - number of points per cell
353  // \endcode
354 
355  // \param baryValues [out] - Output array of barycentric coords
356  // \param cartValues [in] - Input array of Cartesian coords
357  // \param vertices [out] - Vertices of each cell.
358 
359  // */
360  // template<class Scalar, class ArrayTypeOut, class ArrayTypeIn1, class ArrayTypeIn2>
361  // static void cartToBaryTriangle( ArrayTypeOut & baryValues ,
362  // const ArrayTypeIn1 & cartValues ,
363  // const ArrayTypeIn2 & vertices );
364 
365 
366  // /* \brief Converts barycentric coordinates to Cartesian coordinates
367  // on a batch of triangles.
368  // The input array baryValues is (C,P,3)
369  // The output array cartValues is (C,P,2).
370  // The input array vertices is (C,3,2), where
371  // \code
372  // C - num. integration domains
373  // P - number of points per cell
374  // D - is the spatial dimension
375  // \endcode
376 
377  // \param baryValues [out] - Output array of barycentric coords
378  // \param cartValues [in] - Input array of Cartesian coords
379  // \param vertices [out] - Vertices of each cell.
380 
381  // */
382  // template<class Scalar, class ArrayTypeOut, class ArrayTypeIn1, class ArrayTypeIn2>
383  // static void baryToCartTriangle( ArrayTypeOut & cartValues ,
384  // const ArrayTypeIn1 & baryValues ,
385  // const ArrayTypeIn2 & vertices );
386 
387 
388  // /* \brief Converts Cartesian coordinates to barycentric coordinates
389  // on a batch of tetrahedra.
390  // The input array cartValues is (C,P,3)
391  // The output array baryValues is (C,P,4).
392  // The input array vertices is (C,4,3), where
393  // \code
394  // C - num. integration domains
395  // P - number of points per cell
396  // D - is the spatial dimension
397  // \endcode
398 
399  // \param baryValues [out] - Output array of barycentric coords
400  // \param cartValues [in] - Input array of Cartesian coords
401  // \param vertices [out] - Vertices of each cell.
402 
403  // */
404  // template<class Scalar, class ArrayTypeOut, class ArrayTypeIn1, class ArrayTypeIn2>
405  // static void cartToBaryTetrahedron( ArrayTypeOut & baryValues ,
406  // const ArrayTypeIn1 & cartValues ,
407  // const ArrayTypeIn2 & vertices );
408 
409  // /* \brief Converts barycentric coordinates to Cartesian coordinates
410  // on a batch of tetrahedra.
411  // The input array baryValues is (C,P,4)
412  // The output array cartValues is (C,P,3).
413  // The input array vertices is (C,4,3), where
414  // \code
415  // C - num. integration domains
416  // P - number of points per cell
417  // D - is the spatial dimension
418  // \endcode
419 
420  // \param baryValues [out] - Output array of barycentric coords
421  // \param cartValues [in] - Input array of Cartesian coords
422  // \param vertices [out] - Vertices of each cell.
423 
424  // */
425  // template<class Scalar, class ArrayTypeOut, class ArrayTypeIn1, class ArrayTypeIn2>
426  // static void baryToCartTetrahedron( ArrayTypeOut & cartValues ,
427  // const ArrayTypeIn1 & baryValues ,
428  // const ArrayTypeIn2 & vertices );
429 
430 
431 
432 
447  template<typename pointValueType, class ...pointProperties>
448  static void
449  getEquispacedLatticeTriangle( Kokkos::DynRankView<pointValueType,pointProperties...> points,
450  const ordinal_type order,
451  const ordinal_type offset = 0 );
452 
466  template<typename pointValueType, class ...pointProperties>
467  static void
468  getEquispacedLatticeTetrahedron( Kokkos::DynRankView<pointValueType,pointProperties...> points,
469  const ordinal_type order ,
470  const ordinal_type offset = 0 );
471 
472 
479  template<typename pointValueType, class ...pointProperties>
480  static void
481  warpFactor( Kokkos::DynRankView<pointValueType,pointProperties...> warp ,
482  const ordinal_type order ,
483  const Kokkos::DynRankView<pointValueType,pointProperties...> xnodes ,
484  const Kokkos::DynRankView<pointValueType,pointProperties...> xout
485  );
486 
501  template<typename pointValueType, class ...pointProperties>
502  static void
503  getWarpBlendLatticeTriangle( Kokkos::DynRankView<pointValueType,pointProperties...> points,
504  const ordinal_type order ,
505  const ordinal_type offset = 0 );
506 
520  template<typename pointValueType, class ...pointProperties>
521  static void
522  getWarpBlendLatticeTetrahedron( Kokkos::DynRankView<pointValueType,pointProperties...> points,
523  const ordinal_type order ,
524  const ordinal_type offset = 0 );
525 
526 
537  template<typename pointValueType, class ...pointProperties>
538  static void
539  warpShiftFace3D( Kokkos::DynRankView<pointValueType,pointProperties...> dxy,
540  const ordinal_type order ,
541  const pointValueType pval ,
542  const Kokkos::DynRankView<pointValueType,pointProperties...> L1,
543  const Kokkos::DynRankView<pointValueType,pointProperties...> L2,
544  const Kokkos::DynRankView<pointValueType,pointProperties...> L3,
545  const Kokkos::DynRankView<pointValueType,pointProperties...> L4
546  );
547 
557  template<typename pointValueType, class ...pointProperties>
558  static void
559  evalshift( Kokkos::DynRankView<pointValueType,pointProperties...> dxy ,
560  const ordinal_type order ,
561  const pointValueType pval ,
562  const Kokkos::DynRankView<pointValueType,pointProperties...> L1 ,
563  const Kokkos::DynRankView<pointValueType,pointProperties...> L2 ,
564  const Kokkos::DynRankView<pointValueType,pointProperties...> L3
565  );
566 
573  template<typename pointValueType, class ...pointProperties>
574  static void
575  evalwarp( Kokkos::DynRankView<pointValueType,pointProperties...> warp ,
576  const ordinal_type order ,
577  const Kokkos::DynRankView<pointValueType,pointProperties...> xnodes ,
578  const Kokkos::DynRankView<pointValueType,pointProperties...> xout );
579 };
580 
581 }
582 
584 
585 #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 getWarpBlendLattice(Kokkos::DynRankView< pointValueType, pointProperties...> points, const shards::CellTopology cellType, const ordinal_type order, const ordinal_type offset=0)
Computes a warped lattice (ala Warburton&#39;s warp-blend points of a given order on a reference simplex ...
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 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 (currently disabled for other ce...
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 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...
Utility class that provides methods for calculating distributions of points on different cells...
static void getEquispacedLattice(Kokkos::DynRankView< pointValueType, pointProperties...> points, const shards::CellTopology cellType, const ordinal_type order, const ordinal_type offset=0)
Computes an equispaced lattice of a given order on a reference simplex (currently disabled for other ...
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 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...