Basix 0.9.0

Home     Installation     Demos     C++ docs     Python docs

basix::lattice Namespace Reference

Lattices of points. More...

Enumerations

enum class  type {
  equispaced = 0 , gll = 1 , chebyshev = 2 , gl = 4 ,
  chebyshev_plus_endpoints = 10 , gl_plus_endpoints = 11
}
 The type of point spacing to be used in a lattice. More...
 
enum class  simplex_method { none = 0 , warp = 1 , isaac = 2 , centroid = 3 }
 The method used to generate points inside simplices. More...
 

Functions

template<std::floating_point T>
std::pair< std::vector< T >, std::array< std::size_t, 2 > > create (cell::type celltype, int n, lattice::type type, bool exterior, lattice::simplex_method simplex_method=lattice::simplex_method::none)
 Create a lattice of points on a reference cell optionally including the outer surface points. More...
 

Detailed Description

Lattices of points.

Enumeration Type Documentation

◆ type

enum basix::lattice::type
strong

The type of point spacing to be used in a lattice.

Note
type::chebyshev_plus_endpoints() and type::gl_plus_endpoints() are only intended for internal use only.
Enumerator
equispaced 

Equally spaced points

gll 

Gauss-Lobatto-Legendre (GLL) points

chebyshev 

Chebyshev points

gl 

Gauss-Legendre (GL) points

chebyshev_plus_endpoints 

Chebyshev points plus the endpoints of the interval

gl_plus_endpoints 

Gauss-Legendre (GL) points plus the endpoints of the interval

◆ simplex_method

The method used to generate points inside simplices.

Enumerator
none 

Used when no method is needed, e.g. when making points on a quadrilateral, or when making equispaced points).

warp 

Warping from Hesthaven and Warburton, Nodal Discontinuous Galerkin Methods, https://doi.org/10.1007/978-0-387-72067-8.

isaac 

Points described in Isaac, Recursive, Parameter-Free, Explicitly Defined Interpolation Nodes for Simplices, https://doi.org/10.1137/20M1321802.

centroid 

Place points at the centroids of the grid created by putting points on the edges, as described in Blyth and Pozrikidis, A Lobatto interpolation grid over the triangle, https://doi.org/10.1093/imamat/hxh077.

Function Documentation

◆ create()

template<std::floating_point T>
std::pair< std::vector< T >, std::array< std::size_t, 2 > > basix::lattice::create ( cell::type  celltype,
int  n,
lattice::type  type,
bool  exterior,
lattice::simplex_method  simplex_method = lattice::simplex_method::none 
)

Create a lattice of points on a reference cell optionally including the outer surface points.

For a given celltype, this creates a set of points on a regular grid which covers the cell, eg for a quadrilateral, with n=2, the points are: [0,0], [0.5,0], [1,0], [0,0.5], [0.5,0.5], [1,0.5], [0,1], [0.5,1], [1,1]. If the parameter exterior is set to false, the points lying on the external boundary are omitted, in this case for a quadrilateral with n == 2, the points are: [0.5, 0.5]. The lattice type can be chosen as type::equispaced or type::gll. The type::gll lattice has points spaced along each edge at the Gauss-Lobatto-Legendre quadrature points. These are the same as type::equispaced when n < 3.

Parameters
celltypeThe cell type.
nSize in each direction. There are n + 1 points along each edge of the cell.
typeA lattice type.
exteriorIf set, includes outer boundaries.
simplex_methodThe method used to generate points on simplices.
Returns
Set of points. Shape is (npoints, tdim) and storage is row-major.