Tools

Version query

modepy.version.VERSION

A tuple like (2013, 1, 1) indicating modepy’s version.

modepy.version.VERSION_TEXT

A textual representation of the modepy version.

Interpolation matrices

modepy.vandermonde(functions: Sequence[Callable[[ndarray], ndarray]], nodes: ndarray) ndarray[source]

Return a (generalized) Vandermonde matrix.

The Vandermonde Matrix is given by \(V_{i,j} := f_j(x_i)\) where functions is the list of \(f_j\) and nodes is the array of \(x_i\), shaped as (d, npts), where d is the number of dimensions and npts is the number of nodes.

modepy.multi_vandermonde(functions: Sequence[Callable[[ndarray], Sequence[ndarray]]], nodes: ndarray) tuple[ndarray, ...][source]

Evaluate multiple (generalized) Vandermonde matrices.

The Vandermonde Matrix is given by \(V_{i,j} := f_j(x_i)\) where functions is the list of \(f_j\) and nodes is the array of \(x_i\), shaped as (d, npts), where d is the number of dimensions and npts is the number of nodes. functions must return tuple instances. A sequence of matrices is returned–i.e. this function works directly on modepy.Basis.gradients() and returns a tuple of matrices.

Note

If only one of the matrices is needed, it may be convenient to instead call vandermonde() with the result of Basis.derivatives().

Vandermonde matrices are very useful to concisely express interpolation. For instance, one may use the inverse \(V^{-1}\) of a Vandermonde matrix \(V\) to map nodal values (i.e. point values at the nodes corresponding to \(V\)) into modal (i.e. basis) coefficients. \(V\) itself maps modal coefficients to nodal values. This allows easy resampling:

modepy.resampling_matrix(basis: Sequence[Callable[[ndarray], ndarray]], new_nodes: ndarray, old_nodes: ndarray, least_squares_ok: bool = False) ndarray[source]

Return a matrix that maps nodal values on old_nodes onto nodal values on new_nodes.

Parameters:
  • basis – A sequence of basis functions accepting arrays of shape (dims, npts), like those returned by modepy.orthonormal_basis_for_space().

  • new_nodes – An array of shape (dims, n_new_nodes)

  • old_nodes – An array of shape (dims, n_old_nodes)

  • least_squares_ok – If False, then nodal values at old_nodes are required to determine the interpolant uniquely, i.e. the Vandermonde matrix must be square. If True, then a point-wise least-squares best-approximant is used (by ways of the pseudo-inverse of the Vandermonde matrix).

Vandermonde matrices also obey the following relationship useful for obtaining point interpolants:

\[V^T [\text{interpolation coefficients to point $x$}] = \phi_i(x),\]

where \((\phi_i)_i\) is the basis of functions underlying \(V\).

modepy.inverse_mass_matrix(basis: Basis | Sequence[Callable[[ndarray], ndarray]], nodes: ndarray) ndarray[source]

Return a matrix \(A=M^{-1}\), which is the inverse of the one returned by mass_matrix(). Requires that the basis is orthonormal with weight 1.

Added in version 2015.1.

modepy.mass_matrix(basis: Basis | Sequence[Callable[[ndarray], ndarray]], nodes: ndarray) ndarray[source]

Return a mass matrix \(M\), which obeys

\[M_{ij} = \int_\triangle \phi_i(x) \phi_j(x) dx = (V^{-T} V^{-1})_{ij}.\]
Parameters:

basis – assumed to be an orthonormal basis with respect to the \(L^2\) inner product.

Added in version 2014.1.

modepy.modal_mass_matrix_for_face(face: Face, face_quad: Quadrature, trial_functions: Sequence[Callable[[ndarray], ndarray]], test_functions: Sequence[Callable[[ndarray], ndarray]]) ndarray[source]

Using the quadrature face_quad, provide a matrix \(M_f\) that satisfies:

\[\displaystyle {(M_f)}_{ij} \approx \int_F \phi_i(r) \psi_j(r) dr,\]

where \(\phi_i\) are the (volume) basis functions test_functions, and \(\psi_j\) are the (surface) trial_functions.

Added in version 2020.3.

modepy.nodal_mass_matrix_for_face(face: Face, face_quad: Quadrature, trial_functions: Sequence[Callable[[ndarray], ndarray]], test_functions: Sequence[Callable[[ndarray], ndarray]], volume_nodes: ndarray, face_nodes: ndarray) ndarray[source]

Using the quadrature face_quad, provide a matrix \(M_f\) that satisfies:

\[\displaystyle {(M_f)}_{ij} \approx \int_F \phi_i(r) \psi_j(r) dr,\]

where \(\phi_i\) are the (volume) Lagrange basis functions obtained from test_functions at volume_nodes, and \(\psi_j\) are the (surface) Lagrange basis functions obtained from trial_functions at face_nodes.

Added in version 2020.3.

modepy.modal_quad_mass_matrix_for_face(face: Face, face_quad: Quadrature, test_functions: Sequence[Callable[[ndarray], ndarray]]) ndarray[source]

Using the quadrature face_quad, provide a matrix \(M_f\) that satisfies:

\[\displaystyle (M_f \boldsymbol u)_i = \sum_j w_j \phi_i(r_j) u_j,\]

where \(\phi_i\) are the test_functions, \(w_j\) and \(r_j\) are the weights and nodes from face_quad, and \(u_j\) are the point values of the trial function at the nodes of face_quad.

Added in version 2024.2.

modepy.nodal_quad_mass_matrix_for_face(face: Face, face_quad: Quadrature, test_functions: Sequence[Callable[[ndarray], ndarray]], volume_nodes: ndarray) ndarray[source]

Using the quadrature face_quad, provide a matrix \(M_f\) that satisfies:

\[\displaystyle (M_f \boldsymbol u)_i = \sum_j w_j \phi_i(r_j) u_j,\]

where \(\phi_i\) are the (volume) Lagrange basis functions obtained from test_functions at volume_nodes, \(w_j\) and \(r_j\) are the weights and nodes from face_quad, and \(u_j\) are the point values of the trial function at the nodes of face_quad.

Added in version 2021.2.

modepy.spectral_diag_nodal_mass_matrix(quadrature: TensorProductQuadrature) ndarray[source]

Return the diagonal mass matrix for use in the spectral element method. This mass matrix is exact for Lagrange polynomials with respect to Gauss-Legendre (GL) nodes, using GL nodal degrees of freedom. It is approximate for Lagrange polynomials with respect to the Gauss-Lobatto-Legendre (GLL) nodes, using GLL nodal degrees of freedom.

Returns the vector of diagonal entries.

Added in version 2024.2.

Differentiation is also convenient to express by using \(V^{-1}\) to obtain modal values and then using a Vandermonde matrix for the derivatives of the basis to return to nodal values.

modepy.differentiation_matrices(basis: Sequence[Callable[[ndarray], ndarray]], grad_basis: Sequence[Callable[[ndarray], Sequence[ndarray]]], nodes: ndarray, from_nodes: ndarray | None = None) tuple[ndarray, ...][source]

Return matrices carrying out differentiation on nodal values in the \((r,s,t)\) unit directions. (See Coordinates on the triangle and Coordinates on the tetrahedron.)

Parameters:
  • basis – A sequence of basis functions accepting arrays of shape (dims, npts), like those returned by modepy.Basis.functions().

  • grad_basis – A sequence of functions returning the gradients of basis, like those in modepy.Basis.gradients.

  • nodes – An array of shape (dims, n_nodes)

  • from_nodes – An array of shape (dims, n_from_nodes). If None, assumed to be the same as nodes.

Returns:

If grad_basis returned tuples (i.e. in 2D and 3D), returns a tuple of length dims containing differentiation matrices. If not, returns just one differentiation matrix.

Changed in version 2013.4: Added from_nodes.

modepy.diff_matrices(basis: Basis, nodes: ndarray, from_nodes: ndarray | None = None)[source]

Like differentiation_matrices(), but for a given Basis.

Added in version 2024.2.

modepy.diff_matrix_permutation(node_tuples: Sequence[tuple[int, ...]], ref_axis: int) ndarray[source]

Return a numpy array permutation of integers so that:

diff_matrices[ref_axis] == diff_matrices[0][permutation][:, permutation]

Added in version 2020.1.

Modal decay/residual

Estimate the smoothness of a function represented in a basis returned by modepy.orthonormal_basis_for_space().

The method implemented in this module follows this article:

A. Kloeckner, T. Warburton, and J. S. Hesthaven. “Viscous Shock Capturing in a Time-Explicit Discontinuous Galerkin Method”. Mathematical Modelling of Natural Phenomena 6, No. 3 (May 16, 2011): 57-83. http://dx.doi.org/10.1051/mmnp/20116303 http://arxiv.org/abs/1102.3190

Added in version 2013.2.

modepy.modal_decay.simplex_interp_error_coefficient_estimator_matrix(unit_nodes: ndarray, order: int, n_tail_orders: int) ndarray[source]

Return a matrix \(C\) that, when multiplied by a vector of nodal values, yields the coeffiicients belonging to the basis functions of the n_tail_orders highest orders.

The 2-norm of the resulting coefficients can be used as an estimate of the interpolation error.

Added in version 2018.1.

modepy.modal_decay.fit_modal_decay(coeffs: ndarray, dims: int, n: int, ignored_modes: int = 1) tuple[ndarray, ndarray][source]

Fit a curve to the modal decay on each element.

Parameters:
  • coeffs – an array of shape (num_elements, num_modes) containing modal coefficients of the functions to be analyzed.

  • dims – number of dimensions.

  • ignored_modes – the number of modal coefficients to ignore at the beginning. The default value of ‘1’ ignores the constant.

Returns:

a tuple (expn, constant) of arrays of length num_elements, where the modal decay is fit to the curve constant * total_degree**exponent. -exponent-1 can be used as a rough indicator of how many continuous derivatives the underlying function possesses.

modepy.modal_decay.estimate_relative_expansion_residual(coeffs: ndarray, dims: int, n: int, ignored_modes: int = 1) ndarray[source]

Use the modal fit to estimate the relative residual of the expansion. The arguments to this function exactly match fit_modal_decay().

Returns:

An array of estimates of the fraction of the \(L^2\) norm contained in the (unrepresented) tail of

An idea like this is described in this article:

H. Feng and C. Mavriplis, “Adaptive Spectral Element Simulations of Thin Premixed Flame Sheet Deformations”, Journal of Scientific Computing, Volume 17, Nr. 1, S. 385-395, Dec. 2002. http://dx.doi.org/10.1023/A:1015137722700

For this function, however, the decay curve is fitted using the Kloeckner/Warburton/Hesthaven technique (see above).

Transformations between coordinate systems on the simplex

All of these expect and return arrays of shape (dims, npts).

modepy.tools.equilateral_to_unit(equi)[source]
modepy.tools.barycentric_to_unit(bary)[source]
Parameters:

bary – shaped (dims+1, npoints)

modepy.tools.unit_to_barycentric(unit)[source]
Parameters:

unit – shaped (dims,npoints)

modepy.tools.barycentric_to_equilateral(bary)[source]

Interpolation quality

modepy.tools.estimate_lebesgue_constant(n: int, nodes: ndarray, shape: Shape | None = None, *, visualize: bool = False) float[source]

Estimate the Lebesgue constant of the nodes at polynomial order n.

Parameters:
  • nodes – an array of shape (dim, nnodes) as returned by modepy.warp_and_blend_nodes().

  • shape – a Shape.

  • visualize – visualize the function that gives rise to the returned Lebesgue constant. (2D only for now)

Returns:

the Lebesgue constant, a scalar.

Added in version 2013.2.

Changed in version 2020.2: domain parameter was added with support for nodes on the unit hypercube (i.e. unit square in 2D and unit cube in 3D).

Changed in version 2020.3: Renamed domain to shape.

DOF arrays of tensor product discretizations

modepy.tools.reshape_array_for_tensor_product_space(space: TensorProductSpace, ary: ReshapeableT, axis: int = -1) ReshapeableT[source]

Return a reshaped view of ary that exposes the tensor product nature of the space.

Parameters:
  • ary – an array with axis dimension having a length matching the space_dim of space.

  • axis – an integer that must index a dimension of ary with coefficients corresponding to a tensor-product-structured basis (e.g. modal or nodal coefficients).

Returns:

ary reshaped with axis number axis replaced by a tuple of dimensions matching the dimensions of the spaces making up the tensor product. Variation of the represented function along a given dimension will be represented by variation of array entries along the corresponding array axis.

modepy.tools.unreshape_array_for_tensor_product_space(space: TensorProductSpace, ary: ReshapeableT, axis=-1) ReshapeableT[source]

Undoes the effect of reshape_array_for_tensor_product_space(), given the same space and axis.

Types used in documentation

class modepy.tools.ReshapeableT

An array-like protocol that supports finding the shape and reshaping.