polynomials_on_simplices.polynomial.polynomials_unit_simplex_lagrange_basis module

Polynomials on the m-dimensional unit simplex with values in \(\mathbb{R}^n\), expressed using the Lagrange basis.

\[\begin{split}l(x) = \sum_{\substack{\nu \in \mathbb{N}_0^m \\ |\nu| \leq r}} a_{\nu} l_{\nu, r}(x),\end{split}\]

where \(a_{\nu} \in \mathbb{R}^n\).

Basis polynomials in the Lagrange basis are uniquely determined by selecting a sequence of \(\dim \mathcal{P}_r (\Delta_c^m)\) unique points (Lagrange points) in the unit simplex and demanding that the i:th basis function has the value one at the i:th of these points and zero at all the other points.

Here we have used evenly spaced Lagrange points, so that for \(\dim \mathcal{P}_r (\Delta_c^m)\) we have the Lagrange points

\[\begin{split}\{ x_{\nu} \}_{\substack{\nu \in \mathbb{N}_0^m \\ |\nu| \leq r}}, x_{\nu} = \frac{\nu}{r}.\end{split}\]

The basis polynomials \(l_{\nu, r}(x), \nu \in \mathbb{N}_0^n, |\nu| \leq r\) are thus uniquely determined by the conditions

\[\begin{split}l_{\nu, r}(x_{\mu}) = \begin{cases} 1 & \nu = \mu \\ 0 & \text{else} \end{cases}.\end{split}\]

The set \(\{ l_{\nu, r} \}_{\substack{\nu \in \mathbb{N}_0^m \\ |\nu| \leq r}}\) is a basis for the space of all polynomials of degree less than or equal to r on the unit simplex, \(\mathcal{P}_r (\Delta_c^m)\).

class PolynomialLagrange(coeff, r=None, m=1)[source]

Bases: polynomials_on_simplices.polynomial.polynomials_base.PolynomialBase

Implementation of the abstract polynomial base class for a polynomial on the m-dimensional unit simplex, expressed in the Lagrange basis.

\[l(x) = \sum_{i = 0}^{\dim(\mathcal{P}_r(\mathbb{R}^m)) - 1} a_{\nu_i} l_{\nu_i, r}(x).\]
Parameters:
  • coeff – Coefficients for the polynomial in the Lagrange basis for \(\mathcal{P}_r (\mathbb{R}^m, \mathbb{R}^n). \text{coeff}[i] = a_{\nu_i}\), where \(\nu_i\) is the i:th multi-index in the sequence of all multi-indices of dimension m with norm \(\leq r\) (see polynomials_on_simplices.algebra.multiindex.generate() function). Array of scalars for a scalar valued polynomial (n = 1) and array of n-dimensional vectors for a vector valued polynomial (\(n \geq 2\)).
  • m (int) – Dimension of the domain of the polynomial.
  • r (int) – Degree of the polynomial space. Optional, will be inferred from the number of polynomial coefficients if not specified.
basis()[source]

Get basis for the space \(\mathcal{P}_r (\mathbb{R}^m)\) used to express this polynomial.

Returns:Unique identifier for the basis used.
Return type:str
code_str(fn_name)[source]

Generate a function code string for evaluating this polynomial.

Parameters:fn_name (str) – Name for the function in the generated code.
Returns:Code string for evaluating this polynomial.
Return type:str
degree_elevate(s)[source]

Express the polynomial using a higher degree basis.

Let \(p(x) = \sum_{\substack{\nu \in \mathbb{N}_0^m \\ |\nu| \leq r}} a_{\nu} l_{\nu, r}(x)\) be this polynomial, where \(\{ l_{\nu, r} \}_{\substack{\nu \in \mathbb{N}_0^m \\ |\nu| \leq r}}\) is the Lagrange basis for \(\mathcal{P}_r (\mathbb{R}^m)\). Let \(\{ l_{\nu, s} \}_{\substack{\nu \in \mathbb{N}_0^m \\ |\nu| \leq s}}, s \geq r\) be the Lagrange basis for \(\mathcal{P}_s (\mathbb{R}^m)\). Then this function returns a polynomial \(q(x)\)

\[\begin{split}q(x) = \sum_{\substack{\nu \in \mathbb{N}_0^m \\ |\nu| \leq s}} \tilde{a}_{\nu} l_{\nu, s}(x),\end{split}\]

such that \(p(x) = q(x) \, \forall x \in \Delta_c^m\).

Parameters:s (int) – New degree for the polynomial basis the polynomial should be expressed in.
Returns:Elevation of this polynomial to the higher degree basis.
Return type:PolynomialLagrange.
latex_str()[source]

Generate a Latex string for this polynomial.

Returns:Latex string for this polynomial.
Return type:str
latex_str_expanded()[source]

Generate a Latex string for this polynomial, where each basis function has been expanded in the monomial basis.

Returns:Latex string for this polynomial.
Return type:str
partial_derivative(i=0)[source]

Compute the i:th partial derivative of the polynomial.

Parameters:i (int) – Index of partial derivative.
Returns:i:th partial derivative of this polynomial.
Return type:PolynomialLagrange.
to_monomial_basis()[source]

Compute the monomial representation of this polynomial.

Returns:This polynomial expressed in the monomial basis.
Return type:Polynomial.
dual_lagrange_basis(r, n)[source]

Generate all dual Lagrange base functions for the space \(\mathcal{P}_r(\Delta_c^n)\) (i.e. the Lagrange basis for \(\mathcal{P}_r(\Delta_c^n)^*\)).

See dual_lagrange_basis_fn().

Parameters:
  • n (int) – Dimension of the space.
  • r (int) – Degree of the polynomial space.
Returns:

List of dual base functions.

Return type:

List[callable q(l)].

dual_lagrange_basis_fn(mu, r)[source]

Generate a dual basis function to the Lagrange polynomial basis, i.e. the linear map \(q_{\mu, r} : \mathcal{P}_r(\Delta_c^n) \to \mathbb{R}\) that satisfies

\[q_{\mu, r}(l_{\nu, r}) = \delta_{\mu, \nu},\]

where \(l_{\nu, r}\) is the degree r Lagrange basis polynomial indexed by the multi-index \(\nu\) (see lagrange_basis_fn()) and

\[\begin{split}\delta_{\mu, \nu} = \begin{cases} 1 & \mu = \nu \\ 0 & \text{else} \end{cases}.\end{split}\]
Parameters:
  • mu (int or MultiIndex or Tuple[int, …].) – Multi-index indicating which dual Lagrange basis function should be generated.
  • r (int) – Degree of polynomial space.
Returns:

The dual Lagrange basis function as specified by mu and r.

Return type:

Callable \(q_{\mu, r}(l)\).

dual_vector_valued_lagrange_basis(r, m, n, ordering='interleaved')[source]

Generate all dual Lagrange base functions for the space \(\mathcal{P}_r(\Delta_c^m, \mathbb{R}^n)\) (i.e. the Lagrange basis for \(\mathcal{P}_r(\Delta_c^m, \mathbb{R}^n)^*\)).

See dual_vector_valued_lagrange_basis_fn().

Parameters:
  • m (int) – Dimension of the space.
  • r (int) – Degree of the polynomial space.
  • n (int) – Dimension of the target.
  • ordering (str) – How the vector valued basis functions are ordered. Can be “sequential” or “interleaved”. For sequential, sorting is first done on the index of the component that is non-zero, and then the non-zero component is sorted in the same way as the scalar valued basis functions. For “interleaved” basis functions are first sorted on their non-zero component in the same way as scalar valued basis functions, and then they are sorted on the index of the component that is non-zero.
Returns:

List of dual base functions.

Return type:

List[callable q(l)].

dual_vector_valued_lagrange_basis_fn(mu, r, i, n)[source]

Generate a dual basis function to the vector valued Lagrange polynomial basis, i.e. the linear map \(q_{\mu, r, i} : \mathcal{P}_r(\Delta_c^m, \mathbb{R}^n) \to \mathbb{R}\) that satisfies

\[q_{\mu, r, i}(l_{\nu, r, j}) = \delta_{\mu, \nu} \delta_{i, j},\]

where \(l_{\nu, r, j}\) is the degree r vector valued Lagrange basis polynomial indexed by the multi-index \(\nu\) with a non-zero i:th component (see vector_valued_lagrange_basis_fn()) and

\[\begin{split}\delta_{\mu, \nu} = \begin{cases} 1 & \mu = \nu \\ 0 & \text{else} \end{cases}.\end{split}\]
Parameters:
  • mu (int or MultiIndex or Tuple[int, …].) – Multi-index indicating which dual Lagrange basis function should be generated.
  • r (int) – Degree of polynomial space.
  • i (int) – Integer indicating which dual Lagrange basis function should be generated.
  • n (int) – Dimension of the target.
Returns:

The dual Lagrange basis function as specified by mu, r and i.

Return type:

Callable \(q_{\mu, r, i}(l)\).

generate_lagrange_base_coefficients(r, n)[source]

Generate monomial coefficients for all the Lagrange base polynomials for the space \(\mathcal{P}_r(\Delta_c^n)\), with evenly spaced Lagrange points (see generate_lagrange_points()).

Parameters:
  • n (int) – Dimension of the unit simplex.
  • r (int) – Degree of the polynomial space.
Returns:

Matrix containing the coefficients for each base polynomial as column vectors.

Return type:

Numpy array

generate_lagrange_basis_fn_coefficients(nu, r)[source]

Generate monomial coefficients for a Lagrange basis polynomial in the basis for the space \(\mathcal{P}_r(\Delta_c^n)\), with evenly spaced Lagrange points (see generate_lagrange_points()).

Parameters:
  • nu (int or MultiIndex or Tuple[int, …]) – Multi-index indicating which Lagrange basis polynomial we should generate monomial coefficients for.
  • r (int) – Degree of polynomial.
Returns:

Array containing the coefficients for the basis polynomial.

Return type:

Numpy array

generate_lagrange_point(n, r, nu)[source]

Generate a Lagrange point indexed by a multi-index from the set of evenly spaced Lagrange points on the n-dimensional unit simplex (\(\Delta_c^n\)) (Lagrange basis points are constructed so that each basis function has the value 1 at one of the points, and 0 at all the other points).

\[x_{\nu} = \frac{\nu}{r}.\]
Parameters:
  • n (int) – Dimension of the simplex.
  • r (int) – Degree of the polynomial.
  • nu – Multi-index \(\nu\) indexing the Lagrange point, where \(\frac{\nu_i}{r}\) gives the i:th coordinate of the Lagrange point.
Returns:

Point in the n-dimensional unit simplex.

Return type:

Numpy array

Examples

>>> generate_lagrange_point(1, 2, (1,))
array([0.5])
>>> generate_lagrange_point(2, 2, (1, 0))
array([0.5, 0. ])
>>> generate_lagrange_point(2, 2, (0, 1))
array([0. , 0.5])
generate_lagrange_points(n, r)[source]

Generate evenly spaced Lagrange points on the n-dimensional unit simplex (\(\Delta_c^n\)) (Lagrange basis points are constructed so that each basis function has the value 1 at one of the points, and 0 at all the other points).

\[\begin{split}\{ x_{\nu} \}_{\substack{\nu \in \mathbb{N}_0^n \\ |\nu| \leq r}}, x_{\nu} = \frac{\nu}{r}.\end{split}\]
Parameters:
  • n (int) – Dimension of the simplex.
  • r (int) – Degree of the polynomial.
Returns:

List of points in the n-dimensional unit simplex.

Return type:

Numpy array

Examples

>>> generate_lagrange_points(1, 2)
array([0. , 0.5, 1. ])
>>> generate_lagrange_points(2, 2)
array([[0. , 0. ],
       [0.5, 0. ],
       [1. , 0. ],
       [0. , 0.5],
       [0.5, 0.5],
       [0. , 1. ]])
get_associated_sub_simplex(nu, r, simplex=None)[source]

Get the sub simplex associated with a Lagrange basis polynomial.

For a Lagrange basis polynomial p on a simplex T there exist a unique sub simplex f of T such that

  • \(p|_f \neq 0\),
  • \(p|_g = 0\) for all \(g \in \Delta_k(T), g \neq f\),
  • \(\dim f < \dim h\), where \(h\) is any other sub simplex of T for which the above two conditions hold,

where \(\Delta_k(T)\) is the set of all k-dimensional sub simplices of T and \(k = \dim f\).

Parameters:
  • nu (int or MultiIndex or Tuple[int, …]) – Multi-index indicating for which Lagrange basis polynomial we should get the associated sub simplex.
  • r (int) – Degree of polynomial.
  • simplex (Optional[List[int]]) – Vertex indices for the vertices of the simplex. [0, 1, …, n] assumed if not specified.
Returns:

Tuple containing the associated sub simplex, and the sub multi-index associated with the sub simplex (where all non-zero entries has been removed).

Examples

>>> get_associated_sub_simplex((0,), 1)
([0], (1,))
>>> get_associated_sub_simplex((1,), 1)
([1], (1,))
>>> get_associated_sub_simplex((1, 0, 1), 2, [1, 2, 3, 4])
([2, 4], (1, 1))
get_lagrange_base_coefficients(r, n)[source]

Get monomial coefficients for all the Lagrange base polynomials for the space \(\mathcal{P}_r(\Delta_c^n)\), with evenly spaced Lagrange points (see generate_lagrange_points()).

Parameters:
  • n (int) – Dimension of the unit simplex.
  • r (int) – Degree of the polynomial space.
Returns:

Matrix containing the coefficients for each base polynomial as column vectors.

Return type:

Numpy array

get_lagrange_basis_fn_coefficients(nu, r)[source]

Get monomial coefficients for a Lagrange basis polynomial in the basis for the space \(\mathcal{P}_r(\Delta_c^n)\), with evenly spaced Lagrange points (see generate_lagrange_points()).

Parameters:
  • nu (int or MultiIndex or Tuple[int, …]) – Multi-index indicating which Lagrange basis polynomial we should get monomial coefficients for.
  • r (int) – Degree of polynomial.
Returns:

Array containing the coefficients for the basis polynomial.

Return type:

Numpy array

lagrange_basis(r, n)[source]

Generate all Lagrange base polynomials for the space \(\mathcal{P}_r(\Delta_c^n)\).

Parameters:
  • n (int) – Dimension of the space.
  • r (int) – Degree of the polynomial space.
Returns:

List of base polynomials.

Return type:

List[PolynomialLagrange].

lagrange_basis_fn(nu, r)[source]

Generate a Lagrange basis polynomial on the unit simplex (\(\Delta_c^n\)), where n is equal to the length of nu.

Parameters:
  • nu (int or MultiIndex or Tuple[int, …]) – Multi-index indicating which Lagrange basis polynomial should be generated. The polynomial will have the value 1 at the point associated with the multi-index, and value 0 at all other points.
  • r (int) – Degree of polynomial.
Returns:

The Lagrange base polynomial as specified by nu and r.

Return type:

PolynomialLagrange.

Examples

>>> import sympy as sp
>>> x1, x2 = sp.symbols('x1 x2')
>>> lagrange_basis_fn(1, 1)(x1) - x1
0
>>> sp.simplify(lagrange_basis_fn(2, 2)(x1) - (2*x1**2 - x1))
0
>>> sp.simplify(lagrange_basis_fn((1, 1), 2)((x1, x2)) - 4*x1*x2)
0
lagrange_basis_fn_latex(nu, r)[source]

Generate Latex string for a Lagrange basis polynomial on the unit simplex (\(\Delta_c^n\)), where n is equal to the length of nu.

Parameters:
  • nu (int or MultiIndex or Tuple[int, …]) – Multi-index indicating which Lagrange basis polynomial we should generate Latex string for.
  • r (int) – Degree of polynomial.
Returns:

Latex string for the Lagrange base polynomial as specified by nu and r.

Return type:

str

Examples

>>> lagrange_basis_fn_latex(2, 3)
'-9 / 2 x + 18 x^2 - 27 / 2 x^3'
>>> lagrange_basis_fn_latex((1, 1), 3)
'27 x_1 x_2 - 27 x_1^2 x_2 - 27 x_1 x_2^2'
lagrange_basis_fn_latex_compact(nu, r)[source]

Generate compact Latex string for a Lagrange basis polynomial on the unit simplex (\(\Delta_c^n\)), where n is equal to the length of nu.

Parameters:
  • nu (int or MultiIndex or Tuple[int, …]) – Multi-index indicating which Lagrange basis polynomial we should generate Latex string for.
  • r (int) – Degree of polynomial.
Returns:

Latex string for the Lagrange base polynomial as specified by nu and r.

Return type:

str

Examples

>>> lagrange_basis_fn_latex_compact(2, 3)
'l_{2, 3}(x)'
>>> lagrange_basis_fn_latex_compact((1, 1), 3)
'l_{(1, 1), 3}(x)'
lagrange_basis_fn_monomial(nu, r)[source]

Generate a Lagrange basis polynomial on the unit simplex (\(\Delta_c^n\)), where n is equal to the length of nu, expanded in the monomial basis.

This is the same polynomial as given by the lagrange_basis_fn() function, but expressed in the monomial basis.

Parameters:
  • nu (int or MultiIndex or Tuple[int, …]) – Multi-index indicating which Lagrange basis polynomial should be generated The polynomial will have the value 1 at the point associated with the multi-index, and value 0 at all other points.
  • r (int) – Degree of polynomial.
Returns:

The Lagrange base polynomial as specified by nu and r.

Return type:

Polynomial.

Examples

>>> import sympy as sp
>>> x1, x2 = sp.symbols('x1 x2')
>>> lagrange_basis_fn(1, 1)(x1) - x1
0
>>> sp.simplify(lagrange_basis_fn(2, 2)(x1) - (2*x1**2 - x1))
0
>>> sp.simplify(lagrange_basis_fn((1, 1), 2)((x1, x2)) - 4*x1*x2)
0
lagrange_basis_latex(r, n)[source]

Generate Latex strings for all Lagrange base polynomials for the space \(\mathcal{P}_r(\Delta_c^n)\).

Parameters:
  • n (int) – Dimension of the unit simplex.
  • r (int) – Degree of the polynomial space.
Returns:

List of Latex strings for each Lagrange base polynomial.

Return type:

List[str]

Examples

>>> lagrange_basis_latex(2, 1)
['1 - 3 x + 2 x^2', '4 x - 4 x^2', '-x + 2 x^2']
>>> basis_strings = lagrange_basis_latex(2, 2)
>>> expected_basis_strings = list()
>>> expected_basis_strings.append('1 - 3 x_1 + 2 x_1^2 - 3 x_2 + 4 x_1 x_2 + 2 x_2^2')
>>> expected_basis_strings.append('4 x_1 - 4 x_1^2 - 4 x_1 x_2')
>>> expected_basis_strings.append('-x_1 + 2 x_1^2')
>>> expected_basis_strings.append('4 x_2 - 4 x_1 x_2 - 4 x_2^2')
>>> expected_basis_strings.append('4 x_1 x_2')
>>> expected_basis_strings.append('-x_2 + 2 x_2^2')
>>> basis_strings == expected_basis_strings
True
lagrange_basis_latex_compact(r, n)[source]

Generate compact Latex strings for all Lagrange base polynomials for the space \(\mathcal{P}_r(\Delta_c^n)\).

Parameters:
  • n (int) – Dimension of the unit simplex.
  • r (int) – Degree of the polynomial space.
Returns:

List of Latex strings for each Lagrange base polynomial.

Return type:

List[str]

Examples

>>> lagrange_basis_latex_compact(2, 1)
['l_{0, 2}(x)', 'l_{1, 2}(x)', 'l_{2, 2}(x)']
>>> lagrange_basis_latex_compact(1, 2)
['l_{(0, 0), 1}(x)', 'l_{(1, 0), 1}(x)', 'l_{(0, 1), 1}(x)']
lagrange_basis_monomial(r, n)[source]

Generate all Lagrange base polynomials for the space \(\mathcal{P}_r(\Delta_c^n)\), expanded in the monomial basis.

This is the same set of polynomials as given by the lagrange_basis() function, but expressed in the monomial basis.

Parameters:
  • n (int) – Dimension of the space.
  • r (int) – Degree of the polynomial space.
Returns:

List of base polynomials.

Return type:

List[Polynomial].

unique_identifier_lagrange_basis()[source]

Get unique identifier for the Lagrange polynomial basis on the unit simplex.

Returns:Unique identifier.
Return type:str
unit_polynomial(r=0, m=1, n=1)[source]

Get the Lagrange polynomial \(l \in \mathcal{P}(\Delta_c^m, \mathbb{R}^n)\) which is identically one.

Parameters:
  • m (int) – Dimension of the polynomial domain.
  • n (int) – Dimension of the polynomial target.
  • r (int) – The unit polynomial will be expressed in the Lagrange basis for \(\mathcal{P}_r(\Delta_c^m, \mathbb{R}^n)\).
Returns:

The unit polynomial.

Return type:

PolynomialLagrange.

vector_valued_lagrange_basis(r, m, n, ordering='interleaved')[source]

Generate all Lagrange base polynomials for the space \(\mathcal{P}_r(\Delta_c^m, \mathbb{R}^n)\).

Parameters:
  • m (int) – Dimension of the unit simplex.
  • r (int) – Degree of the polynomial space.
  • n (int) – Dimension of the target.
  • ordering (str) – How the vector valued basis functions are ordered. Can be “sequential” or “interleaved”. For sequential, sorting is first done on the index of the component that is non-zero, and then the non-zero component is sorted in the same way as the scalar valued basis functions. For “interleaved” basis functions are first sorted on their non-zero component in the same way as scalar valued basis functions, and then they are sorted on the index of the component that is non-zero.
Returns:

List of base polynomials.

Return type:

List[PolynomialLagrange].

vector_valued_lagrange_basis_fn(nu, r, i, n)[source]

Generate a vector valued Lagrange basis polynomial on the m-dimensional unit simplex, \(l_{\nu, r, i} : \Delta_c^m \to \mathbb{R}^n\).

The vector valued basis polynomial is generated by specifying a scalar valued basis polynomial and the component of the vector valued basis polynomial that should be equal to the scalar valued basis polynomial. All other components of the vector valued basis polynomial will be zero, i.e.

\[\begin{split}l_{\nu, r, i}^j (x) = \begin{cases} l_{\nu, r} (x), & i = j \\ 0, & \text{else} \end{cases},\end{split}\]

where m is equal to the length of nu.

Parameters:
  • nu (int or MultiIndex or Tuple[int, …]) – Multi-index indicating which scalar valued Lagrange basis polynomial should be generated for the non-zero component.
  • r (int) – Degree of polynomial.
  • i (int) – Index of the vector component that is non-zero.
  • n (int) – Dimension of the target.
Returns:

The Lagrange base polynomial as specified by nu, r, i and n.

Return type:

PolynomialLagrange.

Examples

>>> import sympy as sp
>>> x1, x2 = sp.symbols('x1 x2')
>>> vector_valued_lagrange_basis_fn(0, 1, 0, 2)(x1)
array([-x1 + 1, 0], dtype=object)
>>> vector_valued_lagrange_basis_fn(1, 1, 1, 2)(x1)
array([0, x1], dtype=object)
>>> vector_valued_lagrange_basis_fn((1, 0), 2, 0, 2)((x1, x2))
array([-4*x1**2 - 4*x1*x2 + 4*x1, 0], dtype=object)
>>> vector_valued_lagrange_basis_fn((1, 1), 3, 1, 3)((x1, x2))
array([0, -27*x1**2*x2 - 27*x1*x2**2 + 27*x1*x2, 0], dtype=object)
zero_polynomial(r=0, m=1, n=1)[source]

Get the Lagrange polynomial \(l \in \mathcal{P}(\Delta_c^m, \mathbb{R}^n)\) which is identically zero.

Parameters:
  • m (int) – Dimension of the polynomial domain.
  • n (int) – Dimension of the polynomial target.
  • r (int) – The zero polynomial will be expressed in the Lagrange basis for \(\mathcal{P}_r(\Delta_c^m, \mathbb{R}^n)\).
Returns:

The zero polynomial.

Return type:

PolynomialLagrange.