nodes¶
All of the node sets for the \(d\)-simplex implemented by recursivenodes
.
Node sets on the simplex¶
- recursivenodes.nodes.recursive(d, n, family='lgl', domain='barycentric')[source]¶
Alias of
recursivenodes.recursive_nodes()
.
- recursivenodes.nodes.equispaced(d, n, domain='biunit')[source]¶
Equispaced nodes for polynomials up to degree \(n\) on the \(d\)-simplex.
- Parameters:
d (int) – The dimension of the simplex.
n (int) – The polynomial degree
domain (str, optional) – “The domain optional argument” for the choices and their formal definitions.
- Returns:
Equispaced nodes as a 2D array with one row for each of \(\binom{n+d}{d}\) nodes, and \(d\) columns for coordinates (or \(d+1\) if
domain='barycentric'
).- Return type:
ndarray
Example
>>> import matplotlib.pyplot as plt >>> from numpy import eye >>> from recursivenodes.nodes import equispaced >>> from recursivenodes.utils import coord_map >>> nodes = equispaced(2, 7, domain='equilateral') >>> corners = coord_map(eye(3), 'barycentric', 'equilateral') >>> plt.plot(corners[[0,1,2,0],0], corners[[0,1,2,0],1]) [<matplotlib.lines.Line2D object at ...>] >>> plt.scatter(nodes[:,0], nodes[:,1]) <matplotlib.collections.PathCollection object at ...> >>> plt.gca().set_aspect('equal') >>> plt.title('Equispaced Nodes') Text(0.5, 1.0, 'Equispaced Nodes') >>> plt.show()
- recursivenodes.nodes.equispaced_interior(d, n, domain='biunit')[source]¶
Equispaced nodes for polynomials up to degree \(n\) on the \(d\)-simplex, all of which are in the interior.
- Parameters:
d (int) – The dimension of the simplex.
n (int) – The polynomial degree
domain (str, optional) – The domain of the simplex. See “The domain optional argument” for the choices and their formal definitions.
- Returns:
Equispaced interior nodes as a 2D array with one row for each of \(\binom{n+d}{d}\) nodes, and \(d\) columns for coordinates (or \(d+1\) if
domain='barycentric'
).- Return type:
ndarray
Example
>>> import matplotlib.pyplot as plt >>> from numpy import eye >>> from recursivenodes.nodes import equispaced_interior >>> from recursivenodes.utils import coord_map >>> nodes = equispaced_interior(2, 7, domain='equilateral') >>> corners = coord_map(eye(3), 'barycentric', 'equilateral') >>> plt.plot(corners[[0,1,2,0],0], corners[[0,1,2,0],1]) [<matplotlib.lines.Line2D object at ...>] >>> plt.scatter(nodes[:,0], nodes[:,1]) <matplotlib.collections.PathCollection object at ...> >>> plt.gca().set_aspect('equal') >>> plt.title('Equispaced Interior Nodes') Text(0.5, 1.0, 'Equispaced Interior Nodes') >>> plt.show()
- recursivenodes.nodes.blyth_luo_pozrikidis(d, n, x=None, domain='biunit')[source]¶
Create Blyth-Luo-Pozrikidis nodes from a 1D node set for polynomials up to degree \(n\) on the \(d\)-simplex.
Notes
The Blyth-Luo-Pozrikidis rule places an interior node with multi-index \(\boldsymbol{\alpha}\) at the barycentric point \(\boldsymbol{b}(\boldsymbol{\alpha})\) such that
(1)¶\[b_i(\boldsymbol{\alpha}) = x_{n,\boldsymbol{\alpha}_i} + \frac{1}{d}(1 - \sum_{j\neq i} x_{n,\boldsymbol{\alpha}_j}).\]Points on the boundary look like \((d-1)\)-simplex Blyth-Luo-Pozrikidis nodes.
- Parameters:
d (int) – The dimension of the simplex.
n (int) – The polynomial degree
x (ndarray, optional) – 1D node set on \([0, 1]\) with \(n+1\) points. Lobatto-Gauss-Legendre nodes are used if
x=None
.domain (str, optional) – The domain of the simplex. See “The domain optional argument” for the choices and their formal definitions.
- Returns:
Blyth-Luo-Pozrikidis nodes as a 2D array with one row for each of \(\binom{n+d}{d}\) nodes, and \(d\) columns for coordinates (or \(d+1\) if
domain='barycentric'
).- Return type:
ndarray
Example
This plot shows the Blyth-Luo-Pozrikidis nodes with lines connecting the Lobatto-Gauss-Legendre nodes on the edges. The definition in (1) is designed to place the nodes in the centroids of the triangles created by the intersecting lines.
>>> import matplotlib.pyplot as plt >>> from numpy import eye >>> from recursivenodes.nodes import blyth_luo_pozrikidis >>> from recursivenodes.utils import coord_map, multiindex_equal >>> nodes = blyth_luo_pozrikidis(2, 7, domain='equilateral') >>> corners = coord_map(eye(3), 'barycentric', 'equilateral') >>> # plot grid lines >>> for (i, al) in enumerate(multiindex_equal(3, 7)): ... if min(al) > 0 or max(al) == sum(al): continue ... for (j, be) in enumerate(multiindex_equal(3, 7)): ... if (j <= i or min(be) > 0 or max(be) != max(al) ... or sum(be) != sum(al)): ... continue ... for d in range(3): ... if al[d] == be[d] and al[d] != 0: ... _ = plt.plot(nodes[[i,j],0], nodes[[i,j],1], ... linestyle='--', c='grey') >>> plt.plot(corners[[0,1,2,0],0], corners[[0,1,2,0],1]) [<matplotlib.lines.Line2D object at ...>] >>> plt.scatter(nodes[:,0], nodes[:,1]) <matplotlib.collections.PathCollection object at ...> >>> plt.gca().set_aspect('equal') >>> plt.title('Blyth-Luo-Pozrikidis Nodes') Text(0.5, 1.0, 'Blyth-Luo-Pozrikidis Nodes') >>> plt.show()
References
- recursivenodes.nodes.warburton(d, n, x=None, alpha=None, domain='biunit')[source]¶
Warburton warp & blend nodes from a 1D node set for polynomials up to degree \(n\) on the \(d\)-simplex.
Notes
The Warburton warp & blend nodes define a node’s coordinates as a displacement from equispaced coordinates by blending together distortion maps of the \(d\)-simplex that warp the edge nodes to match the 1D node set. One optimization parameter \(\boldsymbol{\alpha}\) controls the blending in the interior of the simplex.
- Parameters:
d (int) – The dimension of the simplex.
n (int) – The polynomial degree
x (ndarray, optional) – 1D node set on \([0, 1]\) with \(n+1\) points. Lobatto-Gauss-Legendre nodes are used if
x=None
.alpha (float, optional) – The blending parameter. If
alpha=None
, a precomputed optimal parameter is used if known.domain (str, optional) – The domain of the simplex. See “The domain optional argument” for the choices and their formal definitions.
- Returns:
Warburton warp & blend nodes as a 2D array with one row for each of \(\binom{n+d}{d}\) nodes, and \(d\) columns for coordinates (or \(d+1\) if
domain='barycentric'
).- Return type:
ndarray
Example
>>> import matplotlib.pyplot as plt >>> from numpy import eye >>> from recursivenodes.nodes import warburton >>> from recursivenodes.utils import coord_map >>> nodes = warburton(2, 7, domain='equilateral') >>> corners = coord_map(eye(3), 'barycentric', 'equilateral') >>> plt.plot(corners[[0,1,2,0],0], corners[[0,1,2,0],1]) [<matplotlib.lines.Line2D object at ...>] >>> plt.scatter(nodes[:,0], nodes[:,1]) <matplotlib.collections.PathCollection object at ...> >>> plt.gca().set_aspect('equal') >>> plt.title('Warburton Warp & Blend Nodes') Text(0.5, 1.0, 'Warburton Warp & Blend Nodes') >>> plt.show()
References
[War06]
- recursivenodes.nodes.warburton_alpha¶
A dictionary of optimal values of the blending parameter \(\alpha\) computed in [War06]. Keyed by
(d, n)
tuples.Example
>>> warburton_alpha[(2,7)] 1.0999
- recursivenodes.nodes.rapetti_sommariva_vianello(n, domain='unit')[source]¶
Create Rapetti-Sommariva-Vianello nodes for the triangle. that minimize the Lebesgue constant while remaining symmetric and having Lobatto-Gauss-Legendre nodes on the edges.
Notes
The Rapetti-Sommariva-Vianello nodes are implicitly defined by minimizing the Lebesgue constanat while maintaining \(S_3\) symmetry and edge nodes that are Lobatto-Gauss-Legendre nodes.
- Parameters:
n (int) – The polynomial degree, must be \(\leq 18\).
domain (str, optional) – The domain of the simplex. See “The domain optional argument” for the choices and their formal definitions.
- Returns:
Rapetti-Sommariva-Vianellow nodes as a 2D array with one row for each of \(\binom{n+2}{2}\) nodes, and 2 columns for coordinates (or 3 if
domain='barycentric'
).- Return type:
ndarray
Example
This plot shows the Rapetti-Sommariva-Vianello nodes for \(n = 10\), which do not have the lattice-like structure of the explicitly defined nodes in the module.
>>> import matplotlib.pyplot as plt >>> from numpy import eye >>> from recursivenodes.nodes import rapetti_sommariva_vianello >>> from recursivenodes.utils import coord_map >>> nodes = rapetti_sommariva_vianello(10, domain='equilateral') >>> corners = coord_map(eye(3), 'barycentric', 'equilateral') >>> plt.plot(corners[[0,1,2,0],0], corners[[0,1,2,0],1]) [<matplotlib.lines.Line2D object at ...>] >>> plt.scatter(nodes[:,0], nodes[:,1]) <matplotlib.collections.PathCollection object at ...> >>> plt.gca().set_aspect('equal') >>> plt.title('Rapetti-Sommariva-Vianello Nodes') Text(0.5, 1.0, 'Rapetti-Sommariva-Vianello Nodes') >>> plt.show()
References
[RSV12]
Warning
These nodes cannot be used with
recursivenodes.lebesgue.lebesguemax()
for \(n \geq 10\), because they do not have lattice-like structure.
The family
optional argument¶
A full definition of 1D node families is given in the documentation of
recursivenodes.recursive_nodes()
.
Here are ways to specify 1D node families in recursivenodes
:
Shifted Lobatto-Gauss-Jacobi nodes for different values of the weight exponents \(\boldsymbol{\alpha}\) and \(\beta\), which always include the endpoints \(0\) and \(1\):
'lgl'
(str) – Shifted L-G-Legendre (\(\boldsymbol{\alpha} = \beta = 0\)).'lgc'
(str) – Shifted L-G-Chebyshev (\(\boldsymbol{\alpha} = \beta = -1/2\)).('lgg', a)
((str, float)) – Shifted L-G-Gegenbauer (\(\boldsymbol{\alpha} = \beta = a - 1/2\)).
Shifted Gauss-Jacobi nodes for different values of the weight exponents \(\boldsymbol{\alpha}\) and \(\beta\), which are in the interior away from the endpoints \(0\) and \(1\):
'gl'
(str) – Shifted G-Legendre (\(\boldsymbol{\alpha} = \beta = 0\)).'gc'
(str) – Shifted G-Chebyshev (\(\boldsymbol{\alpha} = \beta = -1/2\)).('gg', a)
((str, float)) – Shifted G-Gegenbauer (\(\boldsymbol{\alpha} = \beta = a - 1/2\)).
Equispaced nodes:
'equi'
(str) – Equispaced nodes including the endpoints \(0\) and \(1\).'equi_interior'
(str) – Equispaced nodes excluding the endpoints \(0\) and \(1\). The \(k\)-point equispaced interior nodes are a half-phase off from the \((k+1)\)-point equispaced nodes.>>> equispaced(1, 4, domain='unit') array([[0. ], [0.25], [0.5 ], [0.75], [1. ]]) >>> equispaced_interior(1, 3, domain='unit') array([[0.125], [0.375], [0.625], [0.875]])
The domain
optional argument¶
Four reference domains for the \(d\)-simplex are implemented in recursivenodes
.
'unit'
(str) – In \(\mathbb{R}^d\), defined by \(x_i \geq 0, \sum_i x_i <= 1\).'biunit'
(str) – In \(\mathbb{R}^d\), defined by \(x_i \geq -1, \sum_i x_i <= 2-d\).'equilateral'
(str) – In \(\mathbb{R}^d\), defined by having:edge length 2,
centroid at the origin,
a downward-facing facet,
a projection onto the first \(d-1\) coordinates that is the corresponding equilateral \(d-1\)-simplex.
(For \(d=1\),
'equilateral'
is the same as'biunit'
).'barycentric'
(str) – In \(\mathbb{R}^{d+1}\), defined by \(x_i \geq 0, \sum_i x_i = 1`\).
See also