Intrepid
Intrepid_PointTools.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ************************************************************************
3 //
4 // Intrepid 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 Pavel Bochev (pbboche@sandia.gov)
38 // Denis Ridzal (dridzal@sandia.gov), or
39 // Kara Peterson (kjpeter@sandia.gov)
40 //
41 // ************************************************************************
42 // @HEADER
43 
51 #ifndef INTREPID_POINTTOOLS_HPP
52 #define INTREPID_POINTTOOLS_HPP
53 
54 #include "Shards_CellTopology.hpp"
55 #include "Teuchos_Assert.hpp"
56 #include "Intrepid_Polylib.hpp"
58 #include "Intrepid_CellTools.hpp"
59 #include <stdexcept>
60 
61 namespace Intrepid {
62 
195  class PointTools {
196  public:
214  static inline int getLatticeSize( const shards::CellTopology& cellType ,
215  const int order ,
216  const int offset = 0 )
217  {
218  switch( cellType.getKey() ) {
219  case shards::Tetrahedron<4>::key:
220  case shards::Tetrahedron<8>::key:
221  case shards::Tetrahedron<10>::key:
222  {
223  const int effectiveOrder = order - 4 * offset;
224  if (effectiveOrder < 0) return 0;
225  else return (effectiveOrder+1)*(effectiveOrder+2)*(effectiveOrder+3)/6;
226  }
227  break;
228  case shards::Triangle<3>::key:
229  case shards::Triangle<4>::key:
230  case shards::Triangle<6>::key:
231  {
232  const int effectiveOrder = order - 3 * offset;
233  if (effectiveOrder < 0) return 0;
234  else return (effectiveOrder+1)*(effectiveOrder+2)/2;
235  }
236  break;
237  case shards::Line<2>::key:
238  case shards::Line<3>::key:
239  {
240  const int effectiveOrder = order - 2 * offset;
241  if (effectiveOrder < 0) return 0;
242  else return (effectiveOrder+1);
243  }
244  break;
245  default:
246  TEUCHOS_TEST_FOR_EXCEPTION( true , std::invalid_argument ,
247  ">>> ERROR (Intrepid::PointTools::getLatticeSize): Illegal cell type" );
248  }
249  }
250 
251 
268  template<class Scalar, class ArrayType>
269  static void getLattice( ArrayType &pts ,
270  const shards::CellTopology& cellType ,
271  const int order ,
272  const int offset = 0 ,
273  const EPointType pointType = POINTTYPE_EQUISPACED );
274 
281  template<class Scalar, class ArrayType>
282  static void getGaussPoints( ArrayType &pts ,
283  const int order );
284 
285 
286  private:
302  template<class Scalar, class ArrayTypeOut, class ArrayTypeIn1, class ArrayTypeIn2>
303  static void cartToBaryTriangle( ArrayTypeOut & baryValues ,
304  const ArrayTypeIn1 & cartValues ,
305  const ArrayTypeIn2 & vertices );
306 
307 
324  template<class Scalar, class ArrayTypeOut, class ArrayTypeIn1, class ArrayTypeIn2>
325  static void baryToCartTriangle( ArrayTypeOut & cartValues ,
326  const ArrayTypeIn1 & baryValues ,
327  const ArrayTypeIn2 & vertices );
328 
329 
346  template<class Scalar, class ArrayTypeOut, class ArrayTypeIn1, class ArrayTypeIn2>
347  static void cartToBaryTetrahedron( ArrayTypeOut & baryValues ,
348  const ArrayTypeIn1 & cartValues ,
349  const ArrayTypeIn2 & vertices );
350 
367  template<class Scalar, class ArrayTypeOut, class ArrayTypeIn1, class ArrayTypeIn2>
368  static void baryToCartTetrahedron( ArrayTypeOut & cartValues ,
369  const ArrayTypeIn1 & baryValues ,
370  const ArrayTypeIn2 & vertices );
371 
372 
388  template<class Scalar, class ArrayType>
389  static void getEquispacedLattice( const shards::CellTopology& cellType ,
390  ArrayType &points ,
391  const int order ,
392  const int offset = 0 );
393 
394 
410  template<class Scalar, class ArrayType>
411  static void getWarpBlendLattice( const shards::CellTopology& cellType ,
412  ArrayType &points ,
413  const int order ,
414  const int offset = 0);
415 
416 
431  template<class Scalar, class ArrayType>
432  static void getEquispacedLatticeLine( ArrayType &points ,
433  const int order ,
434  const int offset = 0 );
435 
450  template<class Scalar, class ArrayType>
451  static void getEquispacedLatticeTriangle( ArrayType &points ,
452  const int order ,
453  const int offset = 0 );
454 
469  template<class Scalar, class ArrayType>
470  static void getEquispacedLatticeTetrahedron( ArrayType &points ,
471  const int order ,
472  const int offset = 0 );
473 
488  template<class Scalar, class ArrayType>
489  static void getWarpBlendLatticeLine( ArrayType &points ,
490  const int order ,
491  const int offset = 0 );
492 
499  template<class Scalar, class ArrayType>
500  static void warpFactor( const int order ,
501  const ArrayType &xnodes ,
502  const ArrayType &xout ,
503  ArrayType &warp );
504 
519  template<class Scalar, class ArrayType>
520  static void getWarpBlendLatticeTriangle(ArrayType &points ,
521  const int order ,
522  const int offset = 0 );
523 
537  template<class Scalar, class ArrayType>
538  static void getWarpBlendLatticeTetrahedron( ArrayType &points ,
539  const int order ,
540  const int offset = 0 );
541 
542 
553  template<class Scalar, class ArrayType>
554  static void warpShiftFace3D( const int order ,
555  const Scalar pval ,
556  const ArrayType &L1,
557  const ArrayType &L2,
558  const ArrayType &L3,
559  const ArrayType &L4,
560  ArrayType &dxy);
561 
571  template<class Scalar, class ArrayType>
572  static void evalshift( const int order ,
573  const Scalar pval ,
574  const ArrayType &L1 ,
575  const ArrayType &L2 ,
576  const ArrayType &L3 ,
577  ArrayType &dxy );
578 
585  template<class Scalar, class ArrayType>
586  static void evalwarp( ArrayType &warp ,
587  const int order ,
588  const ArrayType &xnodes ,
589  const ArrayType &xout );
590 
591 
592  }; // end class PointTools
593 
594 } // end namespace Intrepid
595 
596 // include templated definitions
598 
599 #endif
600 
601 #if defined(Intrepid_SHOW_DEPRECATED_WARNINGS)
602 #ifdef __GNUC__
603 #warning "The Intrepid package is deprecated"
604 #endif
605 #endif
606 
static void getWarpBlendLatticeLine(ArrayType &points, const int order, const int offset=0)
Returns the Gauss-Lobatto points of a given order on the reference line [-1,1]. The output array is (...
static void baryToCartTetrahedron(ArrayTypeOut &cartValues, const ArrayTypeIn1 &baryValues, const ArrayTypeIn2 &vertices)
Converts barycentric coordinates to Cartesian coordinates on a batch of tetrahedra. The input array baryValues is (C,P,4) The output array cartValues is (C,P,3). The input array vertices is (C,4,3), where.
static void getGaussPoints(ArrayType &pts, const int order)
Header file for the Intrepid::CellTools class.
static void warpShiftFace3D(const int order, const Scalar pval, const ArrayType &L1, const ArrayType &L2, const ArrayType &L3, const ArrayType &L4, ArrayType &dxy)
This is used internally to compute the tetrahedral warp-blend points one each face.
Header file for utility class to provide multidimensional containers.
static void warpFactor(const int order, const ArrayType &xnodes, const ArrayType &xout, ArrayType &warp)
interpolates Warburton&#39;s warp function on the line
static void evalwarp(ArrayType &warp, const int order, const ArrayType &xnodes, const ArrayType &xout)
Used internally to compute the warp on edges of a triangle in warp-blend points.
static void getWarpBlendLatticeTetrahedron(ArrayType &points, const int order, const int offset=0)
Returns Warburton&#39;s warp-blend points of a given order on the reference tetrahedron. The output array is (P,3), where.
static void cartToBaryTriangle(ArrayTypeOut &baryValues, const ArrayTypeIn1 &cartValues, const ArrayTypeIn2 &vertices)
Converts Cartesian coordinates to barycentric coordinates on a batch of triangles. The input array cartValues is (C,P,2) The output array baryValues is (C,P,3). The input array vertices is (C,3,2), where.
static void evalshift(const int order, const Scalar pval, const ArrayType &L1, const ArrayType &L2, const ArrayType &L3, ArrayType &dxy)
Used internally to evaluate the point shift for warp-blend points on faces of tets.
static void getEquispacedLatticeTetrahedron(ArrayType &points, const int order, const int offset=0)
Computes an equispaced lattice of a given order on the reference tetrahedron. The output array is (P...
static void getEquispacedLatticeTriangle(ArrayType &points, const int order, const int offset=0)
Computes an equispaced lattice of a given order on the reference triangle. The output array is (P...
static void baryToCartTriangle(ArrayTypeOut &cartValues, const ArrayTypeIn1 &baryValues, const ArrayTypeIn2 &vertices)
Converts barycentric coordinates to Cartesian coordinates on a batch of triangles. The input array baryValues is (C,P,3) The output array cartValues is (C,P,2). The input array vertices is (C,3,2), where.
static void getLattice(ArrayType &pts, const shards::CellTopology &cellType, const int order, const int 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 getWarpBlendLattice(const shards::CellTopology &cellType, ArrayType &points, const int order, const int offset=0)
Computes a warped lattice (ala Warburton&#39;s warp-blend points of a given order on a reference simplex ...
static void cartToBaryTetrahedron(ArrayTypeOut &baryValues, const ArrayTypeIn1 &cartValues, const ArrayTypeIn2 &vertices)
Converts Cartesian coordinates to barycentric coordinates on a batch of tetrahedra. The input array cartValues is (C,P,3) The output array baryValues is (C,P,4). The input array vertices is (C,4,3), where.
static void getWarpBlendLatticeTriangle(ArrayType &points, const int order, const int offset=0)
Returns Warburton&#39;s warp-blend points of a given order on the reference triangle. The output array is...
Definition file for utilities for barycentric coordinates and lattices.
static void getEquispacedLattice(const shards::CellTopology &cellType, ArrayType &points, const int order, const int offset=0)
Computes an equispaced lattice of a given order on a reference simplex (currently disabled for other ...
static void getEquispacedLatticeLine(ArrayType &points, const int order, const int offset=0)
Computes an equispaced lattice of a given order on the reference line [-1,1]. The output array is (P...
static int getLatticeSize(const shards::CellTopology &cellType, const int order, const int offset=0)
Computes the number of points in a lattice of a given order on a simplex (currently disabled for othe...
Utility class that provides methods for calculating distributions of points on different cells...
Header file for a set of functions providing orthogonal polynomial polynomial calculus and interpolat...