Connections: Interpolating between Discretizations¶

Base classes¶

class meshmode.discretization.connection.DiscretizationConnection(from_discr: Discretization, to_discr: Discretization, is_surjective: bool)[source]¶

Abstract interface for transporting a DOF vector from one meshmode.discretization.Discretization to another. Possible applications include:

  • upsampling/downsampling on the same mesh

  • restricition to the boundary

  • interpolation to a refined/coarsened mesh

  • interpolation onto opposing faces

  • computing modal data from nodal coefficients

  • computing nodal coefficients from modal data

from_discr¶
to_discr¶
is_surjective¶

A bool indicating whether every output degree of freedom is set by the connection.

abstract __call__(ary: ArrayOrContainerT) ArrayOrContainerT[source]¶

Apply the connection. If applicable, may return a view of the data instead of a copy, i.e. changes to ary may or may not appear in the result returned by this method, and vice versa.

class meshmode.discretization.connection.IdentityDiscretizationConnection(discr: Discretization)[source]¶

A no-op connection from a Discretization to the same discretization that returns the same data unmodified.

class meshmode.discretization.connection.ChainedDiscretizationConnection(connections, from_discr=None)[source]¶

Aggregates multiple DiscretizationConnection instances into a single one.

connections¶
class meshmode.discretization.connection.L2ProjectionInverseDiscretizationConnection(connections, is_surjective=False)[source]¶

Creates an inverse DiscretizationConnection from an existing connection to allow transporting from the original connection’s to_discr to from_discr.

from_discr¶
to_discr¶
is_surjective¶
conn¶
__call__(ary)[source]¶
Parameters:

ary – a DOFArray, or an arraycontext.ArrayContainer of them, containing nodal coefficient data on from_discr.

class meshmode.discretization.connection.DirectDiscretizationConnection(from_discr: Discretization, to_discr: Discretization, groups: Sequence[DiscretizationConnectionElementGroup], is_surjective: bool)[source]¶

A concrete DiscretizationConnection supported by interpolation data.

from_discr¶
to_discr¶
groups¶

a list of DiscretizationConnectionElementGroup instances, with a one-to-one correspondence to the groups in to_discr.

is_surjective¶

A bool indicating whether every output degree of freedom is set by the connection.

is_permutation(tol_multiplier: float | None = None) bool[source]¶

Return True if no interpolation is used in applying this connection, i.e. if the source unit nodes in the connection (cf. InterpolationBatch.result_unit_nodes) match up with the unit nodes in the source element, up to machine epsilon multiplied by tol_multiplier (or an internally specified tolerance if none is given).

__call__(ary: ArrayOrContainerT, *, _force_use_loopy: bool = False, _force_no_merged_batches: bool = False) ArrayOrContainerT[source]¶
Parameters:

ary – a DOFArray, or an arraycontext.ArrayContainer of them, containing nodal coefficient data on from_discr.

Mapping between modal and nodal representations¶

class meshmode.discretization.connection.NodalToModalDiscretizationConnection(from_discr, to_discr, allow_approximate_quad=False)[source]¶

A concrete subclass of DiscretizationConnection, which maps nodal data to its modal representation. For interpolatory (unisolvent) element groups, the mapping from nodal to modal representations is performed via:

\[y = V^{-1} [\text{nodal basis coefficients}]\]

where \(V_{i,j} = \phi_j(x_i)\) is the generalized Vandermonde matrix, \(\phi_j\) is a nodal basis, and \(x_i\) are nodal points on the reference element defining the nodal discretization.

For non-interpolatory element groups (for example, QuadratureSimplexElementGroup), modal coefficients are computed using the underlying quadrature rule \((w_q, x_q)\), and an orthonormal basis \(\psi_i\) spanning the modal discretization space. The modal coefficients are then obtained via:

\[y = V^T W [\text{nodal basis coefficients}]\]

where \(V_{i, j} = \psi_j(x_i)\) is the Vandermonde matrix constructed from the orthonormal basis evaluated at the quadrature nodes \(x_i\), and \(W = \text{Diag}(w_q)\) is a diagonal matrix containing the quadrature weights \(w_q\).

Note

This connection requires that both nodal and modal discretizations are defined on the same mesh.

from_discr¶

An instance of meshmode.discretization.Discretization containing NodalElementGroupBase element groups

to_discr¶

An instance of meshmode.discretization.Discretization containing ModalElementGroupBase element groups

groups¶

A list of DiscretizationConnectionElementGroup instances, with a one-to-one correspondence to the groups in to_discr.

__init__(from_discr, to_discr, allow_approximate_quad=False)[source]¶
Parameters:
__call__(ary)[source]¶

Computes modal coefficients data from a functions nodal coefficients.

Parameters:

ary – a DOFArray, or an arraycontext.ArrayContainer of them, containing nodal coefficient data.

class meshmode.discretization.connection.ModalToNodalDiscretizationConnection(from_discr, to_discr)[source]¶

A concrete subclass of DiscretizationConnection, which maps modal data back to its nodal representation. This is computed via:

\[y = V [\text{modal basis coefficients}]\]

where \(V_{i,j} = \phi_j(x_i)\) is the generalized Vandermonde matrix, \(\phi_j\) is an orthonormal (modal) basis, and \(x_i\) are nodal points on the reference element defining the nodal discretization.

Note

This connection requires that both nodal and modal discretizations are defined on the same mesh.

from_discr¶

An instance of meshmode.discretization.Discretization containing ModalElementGroupBase element groups

to_discr¶

An instance of meshmode.discretization.Discretization containing NodalElementGroupBase element groups

groups¶

A list of DiscretizationConnectionElementGroup instances, with a one-to-one correspondence to the groups in to_discr.

__init__(from_discr, to_discr)[source]¶
Parameters:
__call__(ary)[source]¶

Computes nodal coefficients from modal data.

Parameters:

ary – a DOFArray, or an arraycontext.ArrayContainer of them, containing modal coefficient data.

Same-mesh connections¶

meshmode.discretization.connection.make_same_mesh_connection(actx, to_discr, from_discr)[source]¶

Restriction to faces¶

meshmode.discretization.connection.FACE_RESTR_INTERIOR = <class 'meshmode.discretization.connection.face.FACE_RESTR_INTERIOR'>[source]¶

A special value to pass to meshmode.discretization.connection.make_face_restriction() to produce a discretization consisting of all interior faces in a discretization.

meshmode.discretization.connection.FACE_RESTR_ALL = <class 'meshmode.discretization.connection.face.FACE_RESTR_ALL'>[source]¶

A special value to pass to meshmode.discretization.connection.make_face_restriction() to produce a discretization consisting of all faces (interior and boundary) faces in a discretization.

meshmode.discretization.connection.make_face_restriction(actx, discr, group_factory, boundary_tag, per_face_groups=False)[source]¶

Create a mesh, a discretization and a connection to restrict a function on discr to its values on the edges of element faces denoted by boundary_tag.

Parameters:
  • boundary_tag – The boundary tag for which to create a face restriction. May be FACE_RESTR_INTERIOR to indicate interior faces, or FACE_RESTR_ALL to make a discretization consisting of all (interior and boundary) faces.

  • per_face_groups –

    If True, the resulting discretization is guaranteed to have groups organized as:

    (grp0, face0), (grp0, face1), ... (grp0, faceN),
    (grp1, face0), (grp1, face1), ... (grp1, faceN), ...
    

    each with the elements in the same order as the originating group. If False, volume and boundary groups correspond with each other one-to-one, and an interpolation batch is created per face.

Returns:

a meshmode.discretization.connection.DirectDiscretizationConnection representing the new connection. The new boundary discretization can be obtained from the meshmode.discretization.connection.DirectDiscretizationConnection.to_discr attribute of the return value, and the corresponding new boundary mesh from that.

meshmode.discretization.connection.make_face_to_all_faces_embedding(actx, faces_connection, all_faces_discr, from_discr=None)[source]¶

Return a meshmode.discretization.connection.DiscretizationConnection connecting a discretization containing some faces of a discretization to one containing all faces.

Parameters:
meshmode.discretization.connection.make_opposite_face_connection(actx, volume_to_bdry_conn)[source]¶

Given a boundary restriction connection volume_to_bdry_conn, return a DirectDiscretizationConnection that performs data exchange across opposite faces.

Refinement¶

meshmode.discretization.connection.make_refinement_connection(actx, refiner, coarse_discr, group_factory)[source]¶

Return a meshmode.discretization.connection.DiscretizationConnection connecting coarse_discr to a discretization on the fine mesh.

Parameters:

Flattening a ChainedDiscretizationConnection¶

meshmode.discretization.connection.flatten_chained_connection(actx, connection)[source]¶

Collapse a connection into a direct connection.

If the given connection is already a DirectDiscretizationConnection nothing is done. However, if the connection is a ChainedDiscretizationConnection, a new direct connection is constructed that transports from from_discr to to_discr.

The new direct connection will have a number of groups and batches that is, at worse, the product of all the connections in the chain. For example, if we consider a connection between a discretization and a two-level refinement, both levels will have \(n\) groups and \(m + 1\) batches per group, where \(m\) is the number of subdivisions of an element (exact number depends on implementation details in make_refinement_connection()). However, a direct connection from level \(0\) to level \(2\) will have at worst \(n^2\) groups and each group will have \((m + 1)^2\) batches.

Warning

If a large number of connections is chained, the number of groups and batches can become very large.

Parameters:
Returns:

An instance of DirectDiscretizationConnection.

Implementation details¶

class meshmode.discretization.connection.InterpolationBatch(from_group_index: int, from_element_indices: ArrayT, to_element_indices: ArrayT, result_unit_nodes: ndarray, to_element_face: int | None)[source]¶

One interpolation batch captures how a batch of elements within an element group should be an interpolated. Note that while it’s possible that an interpolation batch takes care of interpolating an entire element group from source to target, that’s not necessarily the case. Consider the case of extracting boundary values of a discretization. For, say, a triangle, at least three different interpolation batches are needed to cover boundary edges that fall onto each of the three edges of the unit triangle.

from_group_index¶

An integer indicating from which element group in the from discretization the data should be interpolated.

from_element_indices¶

An array of dtype/shape element_id_t [nelements]. This contains the (group-local) element index (relative to from_group_index from which this “to” element’s data will be interpolated.

to_element_indices¶

An array of dtype/shape element_id_t [nelements]. This contains the (group-local) element index to which this “to” element’s data will be interpolated.

result_unit_nodes¶

A numpy.ndarray of shape (from_group.dim,to_group.nunit_nodes) storing the coordinates of the nodes (in unit coordinates of the from reference element) from which the node locations of this element should be interpolated.

nelements¶
to_element_face¶

int or None. If this interpolation batch targets interpolation to a face, then this number captures the face number (on all elements referenced by from_element_indices to which this batch interpolates. (Since there is a fixed set of “from” unit nodes per batch, one batch will always go to a single face index.)

Note

This attribute is not required. It exists only to carry along metadata from make_face_restriction() to routines that build upon its output, such as make_opposite_face_connection(). If you are not building or consuming face restrictions, it is safe to leave this unset and/or ignore it. This attribute probably belongs in a subclass, but that refactoring hasn’t happened yet. (Sorry!)

class meshmode.discretization.connection.DiscretizationConnectionElementGroup(batches)[source]¶
batches¶

A list of InterpolationBatch instances.