Modes (Basis functions)

Function Spaces

class modepy.FunctionSpace[source]

An opaque object representing a finite-dimensional function space of functions \(\mathbb{R}^n \to \mathbb{R}\).

property order: int

Polynomial degree of the function space, if any.

abstract property spatial_dim: int

The number of spatial dimensions in which the functions in the space operate (\(n\) in the above definition).

abstract property space_dim: int

The number of dimensions of the function space.

class modepy.TensorProductSpace(bases: Tuple[FunctionSpace, ...])[source]

Bases: FunctionSpace

A function space defined as the tensor product of lower dimensional spaces.

To recover the tensor product structure of degree-of-freedom arrays (nodal or modal) associated with this type of space, see reshape_array_for_tensor_product_space().

bases: Tuple[FunctionSpace, ...]

A tuple of the base spaces that take part in the tensor product.

property order: int

Polynomial degree of the functions in the space, if any.

property spatial_dim: int

The number of spatial dimensions in which the functions in the space operate (\(n\) in the above definition).

property space_dim: int

The number of dimensions of the function space.

class modepy.PN(spatial_dim: int, order: int)[source]

Bases: FunctionSpace

The function space of polynomials with total degree \(N\) = order.

\[P^N:=\operatorname{span}\left\{\prod_{i=1}^d x_i^{n_i}:\sum n_i\le N\right\}.\]
property order: int

Total degree of the polynomials in the space.

property spatial_dim: int

The number of spatial dimensions in which the functions in the space operate (\(n\) in the above definition).

property space_dim: int

The number of dimensions of the function space.

class modepy.QN(spatial_dim: int, order: int)[source]

Bases: TensorProductSpace

The function space of polynomials with maximum degree \(N\) = order:

\[Q^N:=\operatorname{span} \left \{\prod_{i=1}^d x_i^{n_i}:\max n_i\le N\right\}.\]
property order

Maximum degree of the polynomials in the space.

modepy.space_for_shape(shape: Shape, order: int | Tuple[int, ...]) FunctionSpace[source]
modepy.space_for_shape(shape: TensorProductShape, order: int | Tuple[int, ...]) TensorProductSpace
modepy.space_for_shape(shape: Simplex, order: int | Tuple[int, ...]) PN

Return an unspecified instance of FunctionSpace suitable for approximation on shape attaining interpolation error of \(O(h^{\text{order}+1})\).

Parameters:

order – an integer interpolation order or a tuple of orders. Taking a tuple of orders is not supported by all function spaces. A notable exception being TensorProductSpace, which allows defining different orders for each base space.

This functionality provides sets of basis functions for the reference elements in modepy.shapes.

class modepy.modes.RealValueT

TypeVar for basis function inputs and outputs.

alias of TypeVar(‘RealValueT’, ~numpy.ndarray, pymbolic.primitives.Expression, float)

Basis Retrieval

exception modepy.BasisNotOrthonormal[source]

Exception raised when a Basis is expected to be orthonormal.

class modepy.Basis[source]
abstract orthonormality_weight() float[source]
Raises:

BasisNotOrthonormal if the basis does not have a weight, i.e. it is not orthogonal. For now, only scalar orthonormality weights are supported. In the future, this may become a symbolic expression involving symbols "rstuvw".

abstract property mode_ids: Tuple[Hashable, ...]

A tuple of of mode (basis function) identifiers, one for each basis function. Mode identifiers should generally be viewed as opaque. They are hashable. For some bases, they are tuples of length matching the number of dimensions and indicating the order along each reference axis.

abstract property functions: Tuple[Callable[[ndarray], ndarray], ...]

A tuple of (callable) basis functions of length matching mode_ids. Each function takes a vector of \((r, s, t)\) reference coordinates (depending on dimension) as input.

abstract property gradients: Tuple[Callable[[ndarray], Tuple[ndarray, ...]], ...]

A tuple of (callable) basis functions of length matching mode_ids. Each function takes a vector of \((r, s, t)\) reference coordinates (depending on dimension) as input. Each function returns a tuple of derivatives, one per reference axis.

modepy.basis_for_space(space: FunctionSpace, shape: Shape) Basis[source]
modepy.basis_for_space(space: PN, shape: Simplex)
modepy.basis_for_space(space: TensorProductSpace, shape: TensorProductShape)
modepy.orthonormal_basis_for_space(space: FunctionSpace, shape: Shape) Basis[source]
modepy.orthonormal_basis_for_space(space: PN, shape: Simplex)
modepy.orthonormal_basis_for_space(space: TensorProductSpace, shape: TensorProductShape)
modepy.monomial_basis_for_space(space: FunctionSpace, shape: Shape) Basis[source]
modepy.monomial_basis_for_space(space: PN, shape: Simplex)
modepy.monomial_basis_for_space(space: TensorProductSpace, shape: TensorProductShape)

Jacobi polynomials

modepy.jacobi(alpha: float, beta: float, n: int, x: RealValueT) RealValueT[source]

Evaluate Jacobi polynomials of type \((\alpha, \beta)\), with \(\alpha, \beta > -1\), and order n at a vector of points x. The points x must lie on the interval \([-1, 1]\).

The polynomials are normalized to be orthonormal with respect to the Jacobi weight \((1 - x)^\alpha (1 + x)^\beta\).

Observe that choosing \(\alpha = \beta = 0\) will yield the Legendre polynomials.

Returns:

a vector of \(P^{(\alpha, \beta)}_n\) evaluated at all x.

modepy.grad_jacobi(alpha: float, beta: float, n: int, x: RealValueT) RealValueT[source]

Evaluate the derivative of jacobi(), with the same meanings and restrictions for all arguments.

Conversion to Symbolic

modepy.symbolicize_function(f: Callable[[RealValueT], RealValueT | Tuple[RealValueT, ...]], dim: int, ref_coord_var_name: str = 'r') RealValueT | Tuple[RealValueT, ...][source]

For a function f (basis or gradient) returned by one of the functions in this module, return a pymbolic expression representing the same function.

Parameters:

dim – the number of dimensions of the reference element on which basis is defined.

Added in version 2020.2.

Tensor product adapter

class modepy.TensorProductBasis(bases: Sequence[Basis], dims_per_basis: Tuple[int, ...] | None = None)[source]

Bases: Basis

Adapts multiple bases into a tensor product basis.

orthonormality_weight() float[source]
Raises:

BasisNotOrthonormal if the basis does not have a weight, i.e. it is not orthogonal. For now, only scalar orthonormality weights are supported. In the future, this may become a symbolic expression involving symbols "rstuvw".

property bases: Sequence[Basis]

A sequence of Basis objects that are being composed into a tensor-product basis in a higher-dimensional space.

property mode_ids: Tuple[Hashable, ...]

A tuple of of mode (basis function) identifiers, one for each basis function. Mode identifiers should generally be viewed as opaque. They are hashable. For some bases, they are tuples of length matching the number of dimensions and indicating the order along each reference axis.

property functions: Tuple[Callable[[ndarray], ndarray], ...]

A tuple of (callable) basis functions of length matching mode_ids. Each function takes a vector of \((r, s, t)\) reference coordinates (depending on dimension) as input.

property gradients: Tuple[Callable[[ndarray], Tuple[ndarray, ...]], ...]

A tuple of (callable) basis functions of length matching mode_ids. Each function takes a vector of \((r, s, t)\) reference coordinates (depending on dimension) as input. Each function returns a tuple of derivatives, one per reference axis.

PKDO basis functions

modepy.modes.pkdo_2d(order: Tuple[int, int], rs: ndarray) ndarray[source]

Evaluate a 2D orthonormal (with weight 1) polynomial on the unit simplex.

Parameters:
  • order – A tuple (i, j) representing the order of the polynomial.

  • rsrs[0], rs[1] are arrays of \((r,s)\) coordinates. (See Coordinates on the triangle)

Returns:

a vector of values of the same length as the rs arrays.

See the following publications:

  • Proriol, Joseph. “Sur une famille de polynomes á deux variables orthogonaux dans un triangle.” CR Acad. Sci. Paris 245 (1957): 2459-2461.

  • Koornwinder, T. “Two-variable analogues of the classical orthogonal polynomials.” Theory and Applications of Special Functions. 1975, pp. 435-495.

  • Dubiner, Moshe. “Spectral Methods on Triangles and Other Domains.” Journal of Scientific Computing 6, no. 4 (December 1, 1991): 345–390. http://dx.doi.org/10.1007/BF01060030

modepy.modes.grad_pkdo_2d(order: Tuple[int, int], rs: ndarray) Tuple[RealValueT, RealValueT][source]

Evaluate the derivatives of pkdo_2d().

Parameters:
  • order – A tuple (i, j) representing the order of the polynomial.

  • rsrs[0], rs[1] are arrays of \((r, s)\) coordinates. (See Coordinates on the triangle)

Returns:

a tuple of vectors (dphi_dr, dphi_ds), each of the same length as the rs arrays.

See the following publications:

  • Proriol, Joseph. “Sur une famille de polynomes á deux variables orthogonaux dans un triangle.” CR Acad. Sci. Paris 245 (1957): 2459-2461.

  • Koornwinder, T. “Two-variable analogues of the classical orthogonal polynomials.” Theory and Applications of Special Functions. 1975, pp. 435-495.

  • Dubiner, Moshe. “Spectral Methods on Triangles and Other Domains.” Journal of Scientific Computing 6, no. 4 (December 1, 1991): 345–390. http://dx.doi.org/10.1007/BF01060030

modepy.modes.pkdo_3d(order: Tuple[int, int, int], rst: ndarray) ndarray[source]

Evaluate a 2D orthonormal (with weight 1) polynomial on the unit simplex.

Parameters:
  • order – A tuple (i, j, k) representing the order of the polynomial.

  • rsrst[0], rst[1], rst[2] are arrays of \((r, s, t)\) coordinates. (See Coordinates on the tetrahedron)

Returns:

a vector of values of the same length as the rst arrays.

See the following publications:

  • Proriol, Joseph. “Sur une famille de polynomes á deux variables orthogonaux dans un triangle.” CR Acad. Sci. Paris 245 (1957): 2459-2461.

  • Koornwinder, T. “Two-variable analogues of the classical orthogonal polynomials.” Theory and Applications of Special Functions. 1975, pp. 435-495.

  • Dubiner, Moshe. “Spectral Methods on Triangles and Other Domains.” Journal of Scientific Computing 6, no. 4 (December 1, 1991): 345–390. http://dx.doi.org/10.1007/BF01060030

modepy.modes.grad_pkdo_3d(order: Tuple[int, int, int], rst: ndarray) Tuple[RealValueT, RealValueT, RealValueT][source]

Evaluate the derivatives of pkdo_3d().

Parameters:
  • order – A tuple (i, j, k) representing the order of the polynomial.

  • rsrs[0], rs[1], rs[2] are arrays of \((r,s,t)\) coordinates. (See Coordinates on the tetrahedron)

Returns:

a tuple of vectors (dphi_dr, dphi_ds, dphi_dt), each of the same length as the rst arrays.

See the following publications:

  • Proriol, Joseph. “Sur une famille de polynomes á deux variables orthogonaux dans un triangle.” CR Acad. Sci. Paris 245 (1957): 2459-2461.

  • Koornwinder, T. “Two-variable analogues of the classical orthogonal polynomials.” Theory and Applications of Special Functions. 1975, pp. 435-495.

  • Dubiner, Moshe. “Spectral Methods on Triangles and Other Domains.” Journal of Scientific Computing 6, no. 4 (December 1, 1991): 345–390. http://dx.doi.org/10.1007/BF01060030

Monomials

modepy.modes.monomial(order: Tuple[int, ...], rst: ndarray) ndarray[source]

Evaluate the monomial of order order at the points rst.

Parameters:
  • order – A tuple (i, j,…) representing the order of the polynomial.

  • rstrst[0], rst[1] are arrays of \((r, s, ...)\) coordinates. (See Coordinates on the triangle)

modepy.modes.grad_monomial(order: Tuple[int, ...], rst: ndarray) Tuple[RealValueT, ...][source]

Evaluate the derivative of the monomial of order order at the points rst.

Parameters:
  • order – A tuple (i, j,…) representing the order of the polynomial.

  • rstrst[0], rst[1] are arrays of \((r, s, ...)\) coordinates. (See Coordinates on the triangle)

Returns:

a tuple of vectors (dphi_dr, dphi_ds, dphi_dt, ….), each of the same length as the rst arrays.

Added in version 2016.1.