Intrepid
|
Utility class that provides methods for calculating distributions of points on different cells. More...
#include <Intrepid_PointTools.hpp>
Static Public Member Functions | |
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 other cell types). If offset == 0, the lattice will include only include the vertex points if order == 1, and will include edge midpoints if order == 2, and so on. In particular, this is the dimension of polynomials of degree "order" on the given simplex. The offset argument is used to indicate that the layer of points on the boundary is omitted (if offset == 1). For greater offsets, more layers are omitteed. More... | |
template<class Scalar , class ArrayType > | |
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 cell types). The output array is (P,D), where. More... | |
template<class Scalar , class ArrayType > | |
static void | getGaussPoints (ArrayType &pts, const int order) |
Static Private Member Functions | |
template<class Scalar , class ArrayTypeOut , class ArrayTypeIn1 , class ArrayTypeIn2 > | |
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. More... | |
template<class Scalar , class ArrayTypeOut , class ArrayTypeIn1 , class ArrayTypeIn2 > | |
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. More... | |
template<class Scalar , class ArrayTypeOut , class ArrayTypeIn1 , class ArrayTypeIn2 > | |
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. More... | |
template<class Scalar , class ArrayTypeOut , class ArrayTypeIn1 , class ArrayTypeIn2 > | |
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. More... | |
template<class Scalar , class ArrayType > | |
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 cell types). The output array is (P,D), where. More... | |
template<class Scalar , class ArrayType > | |
static void | getWarpBlendLattice (const shards::CellTopology &cellType, ArrayType &points, const int order, const int offset=0) |
Computes a warped lattice (ala Warburton's warp-blend points of a given order on a reference simplex (currently disabled for other cell types). The output array is (P,D), where. More... | |
template<class Scalar , class ArrayType > | |
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,1), where. More... | |
template<class Scalar , class ArrayType > | |
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,2), where. More... | |
template<class Scalar , class ArrayType > | |
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,3), where. More... | |
template<class Scalar , class ArrayType > | |
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 (P,1), where. More... | |
template<class Scalar , class ArrayType > | |
static void | warpFactor (const int order, const ArrayType &xnodes, const ArrayType &xout, ArrayType &warp) |
interpolates Warburton's warp function on the line More... | |
template<class Scalar , class ArrayType > | |
static void | getWarpBlendLatticeTriangle (ArrayType &points, const int order, const int offset=0) |
Returns Warburton's warp-blend points of a given order on the reference triangle. The output array is (P,2), where. More... | |
template<class Scalar , class ArrayType > | |
static void | getWarpBlendLatticeTetrahedron (ArrayType &points, const int order, const int offset=0) |
Returns Warburton's warp-blend points of a given order on the reference tetrahedron. The output array is (P,3), where. More... | |
template<class Scalar , class ArrayType > | |
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. More... | |
template<class Scalar , class ArrayType > | |
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. More... | |
template<class Scalar , class ArrayType > | |
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. More... | |
Utility class that provides methods for calculating distributions of points on different cells.
Simplicial lattices in PointTools are sets of points with certain ordering properties. They are used for defining degrees of freedom for higher order finite elements.
Each lattice has an "order". In general, this is the same as the cardinality of the polynomial space of degree "order". In terms of binomial coefficients, this is binomial(order+d,order) for the simplex in d dimensions. On the line, the size is order+1. On the triangle and tetrahedron, there are (order+1)(order+2)/2 and (order+1)(order+2)(order+3)/6, respectively.
The points are ordered lexicographically from low to high in increasing spatial dimension. For example, the line lattice of order 3 looks like:
x--x--x--x
where "x" denotes a point location. These are ordered from left to right, so that the points are labeled:
0--1--2--3
The triangular lattice of order 3 is
x |\ | \ x x | \ | \ x x x | \ | \ x--x--x--x
The ordering starts in the bottom left and increases first from left to right. The ordering is
9 |\ | \ 8 7 | \ | \ 4 5 6 | \ | \ 0--1--2--3
Tetrahedral lattices are similar but difficult to draw with ASCII art.
Each lattice also has an "offset", which indicates a number of layers of points on the bounary taken away. All of the lattices above have a 0 offest. In Intrepid, typically only offset = 0 or 1 will be used. The offset=1 case is used to generate sets of points properly inside a given simplex. These are used, for example, to construct points internal to an edge or face for H(curl) and H(div) finite elements.
For example, for a line lattice with order = 3 and offset = 1, the points will look like
---x--x---
and a triangle with order=3 and offset=1 will contain a single point
. |\ | \ | \ | \ | \ | x \ | \ | \ |--------\
When points on lattices with nonzero offset are numbered, the are numbered contiguously from 0, so that the line and triangle above are respectively
---0--1---
. |\ | \ | \ | \ | \ | 0 \ | \ | \ |--------\
Additionally, two types of point distributions are currently support. The points may be on an equispaced lattice, which is easy to compute but can lead to numerical ill-conditioning in finite element bases and stiffness matrices. Alternatively, the warp-blend points of Warburton are provided on each lattice (which are just the Gauss-Lobatto points on the line).
Definition at line 195 of file Intrepid_PointTools.hpp.
|
staticprivate |
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.
baryValues | [out] - Output array of barycentric coords |
cartValues | [in] - Input array of Cartesian coords |
vertices | [out] - Vertices of each cell. |
|
staticprivate |
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.
baryValues | [out] - Output array of barycentric coords |
cartValues | [in] - Input array of Cartesian coords |
vertices | [out] - Vertices of each cell. |
|
staticprivate |
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.
baryValues | [out] - Output array of barycentric coords |
cartValues | [in] - Input array of Cartesian coords |
vertices | [out] - Vertices of each cell. |
|
staticprivate |
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.
baryValues | [out] - Output array of barycentric coords |
cartValues | [in] - Input array of Cartesian coords |
vertices | [out] - Vertices of each cell. |
|
staticprivate |
Used internally to evaluate the point shift for warp-blend points on faces of tets.
order | [in] - the order of lattice |
pval | [in] - the "alpha" term in the warping function |
L1 | [in] - the first barycentric coordinate of the input points |
L2 | [in] - the second barycentric coordinate of the input points |
L3 | [in] - the third barycentric coordinate of the input points |
dxy | [out] - contains the amount to shift each point in the x and y direction |
Definition at line 491 of file Intrepid_PointToolsDef.hpp.
|
staticprivate |
Used internally to compute the warp on edges of a triangle in warp-blend points.
warp | [out] - a 1d array containing the amount to move each point |
order | [in] - the order of the lattice |
xnodes | [in] - the points to warp to, typically the Gauss-Lobatto points |
xout | [in] - the equispaced points on the edge |
Definition at line 570 of file Intrepid_PointToolsDef.hpp.
References Intrepid::FieldContainer< Scalar, ArrayTypeId >::initialize().
|
staticprivate |
Computes an equispaced lattice of a given order on a reference simplex (currently disabled for other cell types). The output array is (P,D), where.
points | [out] - Output array of point coords |
order | [in] - number of points per side, plus 1 |
offset | [in] - Number of points on boundary to skip |
cellType | [in] - type of reference cell (currently only supports the simplex) |
Definition at line 97 of file Intrepid_PointToolsDef.hpp.
References getLatticeSize().
|
staticprivate |
Computes an equispaced lattice of a given order on the reference line [-1,1]. The output array is (P,1), where.
points | [out] - Output array of point coords |
order | [in] - The lattice has order + 1 points, minus any skipped by offset |
offset | [in] - Number of points on boundary to skip coming in per boundary |
Definition at line 183 of file Intrepid_PointToolsDef.hpp.
|
staticprivate |
Computes an equispaced lattice of a given order on the reference tetrahedron. The output array is (P,3), where.
points | [out] - Output array of point coords |
order | [in] - The lattice has order + 1 points, minus any skipped by offset |
offset | [in] - Number of points on boundary to skip coming in from boundary |
Definition at line 228 of file Intrepid_PointToolsDef.hpp.
|
staticprivate |
Computes an equispaced lattice of a given order on the reference triangle. The output array is (P,2), where.
points | [out] - Output array of point coords |
order | [in] - The lattice has order + 1 points, minus any skipped by offset |
offset | [in] - Number of points on boundary to skip coming in from boundary |
Definition at line 205 of file Intrepid_PointToolsDef.hpp.
|
static |
Retrieves the Gauss-Legendre points from PolyLib, but lets us do it in an arbitrary ArrayType.
pts | [out] - Output array of point coords (P,) |
order | [out] - number of Gauss points - 1 |
Definition at line 79 of file Intrepid_PointToolsDef.hpp.
References Intrepid::IntrepidPolylib::zwgj().
|
static |
Computes a lattice of points of a given order on a reference simplex (currently disabled for other cell types). The output array is (P,D), where.
pts | [out] - Output array of point coords |
cellType | [in] - type of reference cell (currently only supports the simplex) |
order | [in] - number of points per side, plus 1 |
pointType | [in] - flag for point distribution. Currently equispaced and warp/blend points are supported |
offset | [in] - Number of points on boundary to skip |
Definition at line 56 of file Intrepid_PointToolsDef.hpp.
|
inlinestatic |
Computes the number of points in a lattice of a given order on a simplex (currently disabled for other cell types). If offset == 0, the lattice will include only include the vertex points if order == 1, and will include edge midpoints if order == 2, and so on. In particular, this is the dimension of polynomials of degree "order" on the given simplex. The offset argument is used to indicate that the layer of points on the boundary is omitted (if offset == 1). For greater offsets, more layers are omitteed.
cellType | [in] - type of reference cell (currently only supports the simplex) |
order | [in] - order of the lattice |
offset | [in] - the number of boundary layers to omit |
Definition at line 214 of file Intrepid_PointTools.hpp.
Referenced by Intrepid::Basis_HCURL_TET_In_FEM< Scalar, ArrayScalar >::Basis_HCURL_TET_In_FEM(), Intrepid::Basis_HCURL_TRI_In_FEM< Scalar, ArrayScalar >::Basis_HCURL_TRI_In_FEM(), Intrepid::Basis_HDIV_TET_In_FEM< Scalar, ArrayScalar >::Basis_HDIV_TET_In_FEM(), Intrepid::Basis_HDIV_TRI_In_FEM< Scalar, ArrayScalar >::Basis_HDIV_TRI_In_FEM(), getEquispacedLattice(), getWarpBlendLattice(), Intrepid::Basis_HGRAD_TET_Cn_FEM< Scalar, ArrayScalar >::initializeTags(), Intrepid::Basis_HCURL_TET_In_FEM< Scalar, ArrayScalar >::initializeTags(), and main().
|
staticprivate |
Computes a warped lattice (ala Warburton's warp-blend points of a given order on a reference simplex (currently disabled for other cell types). The output array is (P,D), where.
points | [out] - Output array of point coords |
order | [in] - number of points per side, plus 1 |
offset | [in] - Number of points on boundary to skip |
cellType | [in] - type of reference cell (currently only supports the simplex) |
Definition at line 138 of file Intrepid_PointToolsDef.hpp.
References getLatticeSize().
|
staticprivate |
Returns the Gauss-Lobatto points of a given order on the reference line [-1,1]. The output array is (P,1), where.
points | [out] - Output array of point coords |
order | [in] - The lattice has order + 1 points, minus any skipped by offset |
offset | [in] - Number of points on boundary to skip coming in per boundary |
Definition at line 254 of file Intrepid_PointToolsDef.hpp.
References Intrepid::IntrepidPolylib::zwglj().
|
staticprivate |
Returns Warburton's warp-blend points of a given order on the reference tetrahedron. The output array is (P,3), where.
points | [out] - Output array of point coords |
order | [in] - The lattice has order + 1 points, minus any skipped by offset |
offset | [in] - Number of points on boundary to skip coming in per boundary |
Definition at line 616 of file Intrepid_PointToolsDef.hpp.
References Intrepid::FieldContainer< Scalar, ArrayTypeId >::initialize(), Intrepid::CellTools< Scalar >::mapToReferenceFrame(), and Intrepid::FieldContainer< Scalar, ArrayTypeId >::resize().
|
staticprivate |
Returns Warburton's warp-blend points of a given order on the reference triangle. The output array is (P,2), where.
points | [out] - Output array of point coords |
order | [in] - The lattice has order + 1 points, minus any skipped by offset |
offset | [in] - Number of points on boundary to skip coming in per boundary |
Definition at line 336 of file Intrepid_PointToolsDef.hpp.
References Intrepid::CellTools< Scalar >::mapToReferenceFrame().
|
staticprivate |
interpolates Warburton's warp function on the line
order | [in] - The polynomial order |
xnodes | [in] - vector of node locations to interpolate |
xout | [in] - warpfunction at xout, +/- 1 roots deflated |
warp | [out] - the amount to warp each point |
Definition at line 276 of file Intrepid_PointToolsDef.hpp.
References Intrepid::FieldContainer< Scalar, ArrayTypeId >::initialize().
|
staticprivate |
This is used internally to compute the tetrahedral warp-blend points one each face.
order | [in] - the order of lattice |
pval | [in] - the "alpha" term in the warping function |
L1 | [in] - the first barycentric coordinate of the input points |
L2 | [in] - the second barycentric coordinate of the input points |
L3 | [in] - the third barycentric coordinate of the input points |
L4 | [in] - the fourth barycentric coordinate of the input points |
dxy | [out] - contains the amount to shift each point in the x and y direction |
Definition at line 478 of file Intrepid_PointToolsDef.hpp.