Connections: Interpolating between Discretizations¶
Base classes¶
-
class
meshmode.discretization.connection.
DiscretizationConnection
(from_discr, to_discr, is_surjective)¶ 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
-
from_discr
¶
-
to_discr
¶
-
__call__
(queue, vec)¶ Call self as a function.
-
class
meshmode.discretization.connection.
ChainedDiscretizationConnection
(connections, from_discr=None)¶ Aggregates multiple
DiscretizationConnection
instances into a single one.-
connections
¶
-
-
class
meshmode.discretization.connection.
L2ProjectionInverseDiscretizationConnection
(conn, is_surjective=False)¶ 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__
(queue, vec)¶ Call self as a function.
-
-
class
meshmode.discretization.connection.
DirectDiscretizationConnection
(from_discr, to_discr, groups, is_surjective)¶ 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 into_discr
.
-
__call__
(queue, vec)¶ Call self as a function.
-
full_resample_matrix
(queue)¶ Build a dense matrix representing this discretization connection.
Warning
On average, this will be exceedingly expensive (\(O(N^2)\) in the number N of discretization points) in terms of memory usage and thus not what you’d typically want.
-
Same-mesh connections¶
-
meshmode.discretization.connection.
make_same_mesh_connection
(to_discr, from_discr)¶
Restriction to faces¶
-
meshmode.discretization.connection.
FACE_RESTR_INTERIOR
= <class 'meshmode.discretization.connection.face.FACE_RESTR_INTERIOR'>¶ 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'>¶ 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
(discr, group_factory, boundary_tag, per_face_groups=False)¶ 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
meshmode.discretization.connection.FRESTR_INTERIOR_FACES
to indicate interior faces, ormeshmode.discretization.connection.FRESTR_ALL_FACES
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 themeshmode.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
(faces_connection, all_faces_discr, from_discr=None)¶ Return a
meshmode.discretization.connection.DiscretizationConnection
connecting a discretization containing some faces of a discretization to one containing all faces.- Parameters
faces_connection – must be the (connection) result of calling
meshmode.discretization.connection.make_face_restriction()
withmeshmode.discretization.connection.FACE_RESTR_INTERIOR
or a boundary tag.all_faces_discr – must be the (discretization) result of calling
meshmode.discretization.connection.make_face_restriction()
withmeshmode.discretization.connection.FACE_RESTR_ALL
for the same volume discretization as the one from which faces_discr was obtained.from_discr – Allows substituting in a different origin discretization for the returned connection. This discretization must use the same mesh as
faces_connection.to_discr
.
-
meshmode.discretization.connection.
make_opposite_face_connection
(volume_to_bdry_conn)¶ Given a boundary restriction connection volume_to_bdry_conn, return a
DirectDiscretizationConnection
that performs data exchange across opposite faces.
Mesh partitioning¶
-
meshmode.discretization.connection.
make_partition_connection
(local_bdry_conn, i_local_part, remote_bdry, remote_adj_groups, remote_from_elem_faces, remote_from_elem_indices)¶ Connects
local_bdry_conn
to a neighboring partition.- Parameters
local_bdry_conn – A
DiscretizationConnection
of the local partition.i_local_part – The partition number of the local partition.
remote_adj_groups – A list of
InterPartitionAdjacency`
of the remote partition.remote_bdry – A
Discretization
of the boundary of the remote partition.remote_from_elem_faces – remote_from_elem_faces[igrp][idx] gives the face that batch idx interpolates from in group igrp.
remote_from_elem_indices – remote_from_elem_indices[igrp][idx] gives a
np.array
of element indices that batch idx interpolates from in group igrp.
- Returns
A
DirectDiscretizationConnection
that performs data exchange across faces from the remote partition to partition i_local_part.
New in version 2017.1.
Warning
Interface is not final.
Refinement¶
-
meshmode.discretization.connection.
make_refinement_connection
(refiner, coarse_discr, group_factory)¶ Return a
meshmode.discretization.connection.DiscretizationConnection
connecting coarse_discr to a discretization on the fine mesh.- Parameters
refiner – An instance of
meshmode.mesh.refinement.Refiner
coarse_discr – An instance of
meshmode.mesh.discretization.Discretization
associated with the mesh given to the refinergroup_factory – An instance of
meshmode.mesh.discretization.ElementGroupFactory
. Used for discretizing the fine mesh.
Flattening a ChainedDiscretizationConnection
¶
-
meshmode.discretization.connection.
flatten_chained_connection
(queue, connection)¶ Collapse a connection into a direct connection.
If the given connection is already a
DirectDiscretizationConnection
nothing is done. However, if the connection is aChainedDiscretizationConnection
, a new direct connection is constructed that transports fromconnection.from_discr
toconnection.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
queue – An instance of
pyopencl.CommandQueue
.connection – An instance of
DiscretizationConnection
.
- Returns
An instance of
DirectDiscretizationConnection
.
Implementation details¶
-
class
meshmode.discretization.connection.
InterpolationBatch
(from_group_index, from_element_indices, to_element_indices, result_unit_nodes, to_element_face)¶ 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
¶ element_id_t [nelements]
. (apyopencl.array.Array
) This contains the (group-local) element index (relative tofrom_group_index
from which this “to” element’s data will be interpolated.
-
to_element_indices
¶ element_id_t [nelements]
. (apyopencl.array.Array
) 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. (a
pyopencl.array.Array
if existent) If this interpolation batch targets interpolation to a face, then this number captures the face number (on all elements referenced byfrom_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.)
-
-
class
meshmode.discretization.connection.
DiscretizationConnectionElementGroup
(batches)¶ -
batches
¶ A list of
InterpolationBatch
instances.
-