Discretizations

Abstract interface

class meshmode.discretization.ElementGroupBase(mesh_el_group, order, node_nr_base)

Container for the Discretization data corresponding to one meshmode.mesh.MeshElementGroup.

mesh_el_group
order
node_nr_base
nelements
nunit_nodes
nnodes
dim
view(global_array)

Return a view of global_array of shape (..., nelements, nunit_nodes) where global_array is of shape (..., nnodes), where nnodes is the global (per-discretization) node count.

unit_nodes()

Returns a numpy.ndarray of shape (dim, nunit_nodes) of reference coordinates of interpolation nodes.

weights()

Returns an array of length nunit_nodes containing quadrature weights.

class meshmode.discretization.InterpolatoryElementGroupBase(mesh_el_group, order, node_nr_base)

A subclass of ElementGroupBase that is equipped with a function space.

basis()

Returns a list of basis functions that take arrays of shape (dim, n) and return an array of shape (n,)`` (which performs evaluation of the basis function).

grad_basis()
Returns:a tuple of functions, each of which accepts arrays of shape (dims, npts) and returns a tuple of length dims containing the derivatives along each axis as an array of size npts. ‘Scalar’ evaluation, by passing just one vector of length dims, is also supported.
diff_matrices()

Return a dim-long tuple of matrices of shape (nunit_nodes, nunit_nodes), each of which, when applied to an array of nodal values, take derivatives in the reference (r,s,t) directions.

class meshmode.discretization.ElementGroupFactory
__call__(mesh_ele_group, node_nr_base)
class meshmode.discretization.Discretization(cl_ctx, mesh, group_factory, real_dtype=<class 'numpy.float64'>)

An unstructured composite discretization.

Parameters:
real_dtype
complex_dtype
mesh
dim
ambient_dim
nnodes
groups
empty(queue=None, dtype=None, extra_dims=None, allocator=None)

Return an empty DOF vector.

Parameters:dtype – type special value ‘c’ will result in a vector of dtype self.complex_dtype. If None (the default), a real vector will be returned.
zeros(queue, dtype=None, extra_dims=None, allocator=None)
nodes()

shape: (ambient_dim, nnodes)

num_reference_derivative(queue, ref_axes, vec)
quad_weights(queue)

shape: (nnodes)

Composite polynomial discretization

Group types

class meshmode.discretization.poly_element.InterpolatoryQuadratureSimplexElementGroup(mesh_el_group, order, node_nr_base)

Elemental discretization supplying a high-order quadrature rule with a number of nodes matching the number of polynomials in \(P^k\), hence usable for differentiation and interpolation.

No interpolation nodes are present on the boundary of the simplex.

class meshmode.discretization.poly_element.QuadratureSimplexElementGroup(mesh_el_group, order, node_nr_base)

Elemental discretization supplying a high-order quadrature rule with a number of nodes which does not necessarily match the number of polynomials in \(P^k\). This discretization therefore excels at quadarature, but is not necessarily usable for interpolation.

No interpolation nodes are present on the boundary of the simplex.

class meshmode.discretization.poly_element.PolynomialWarpAndBlendElementGroup(mesh_el_group, order, node_nr_base)

Elemental discretization with a number of nodes matching the number of polynomials in \(P^k\), hence usable for differentiation and interpolation. Interpolation nodes edge-clustered for avoidance of Runge phenomena. Nodes are present on the boundary of the simplex.

class meshmode.discretization.poly_element.PolynomialEquidistantSimplexElementGroup(mesh_el_group, order, node_nr_base)

Elemental discretization with a number of nodes matching the number of polynomials in \(P^k\), hence usable for differentiation and interpolation. Interpolation nodes are present on the boundary of the simplex.

New in version 2016.1.

class meshmode.discretization.poly_element.LegendreGaussLobattoTensorProductElementGroup(mesh_el_group, order, node_nr_base)

Group factories

class meshmode.discretization.poly_element.InterpolatoryQuadratureSimplexGroupFactory(order)
class meshmode.discretization.poly_element.QuadratureSimplexGroupFactory(order)
class meshmode.discretization.poly_element.PolynomialWarpAndBlendGroupFactory(order)
class meshmode.discretization.poly_element.PolynomialEquidistantSimplexGroupFactory(order)

New in version 2016.1.

class meshmode.discretization.poly_element.LegendreGaussLobattoTensorProductGroupFactory(order)

Connection/interpolation

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
is_surjective

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

__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.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 in to_discr.

is_surjective

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

__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.

meshmode.discretization.connection.make_same_mesh_connection(to_discr, from_discr)
meshmode.discretization.connection.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()

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, or meshmode.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 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(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:
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.

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_facesremote_from_elem_faces[igrp][idx] gives the face that batch idx interpolates from in group igrp.
  • remote_from_elem_indicesremote_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.

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 refiner
  • group_factory – An instance of meshmode.mesh.discretization.ElementGroupFactory. Used for discretizing the fine mesh.
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 a ChainedDiscretizationConnection, a new direct connection is constructed that transports from connection.from_discr to connection.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, 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]. (a pyopencl.array.Array) 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

element_id_t [nelements]. (a pyopencl.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 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.)

class meshmode.discretization.connection.DiscretizationConnectionElementGroup(batches)
batches

A list of InterpolationBatch instances.

Visualization

meshmode.discretization.visualization.make_visualizer(queue, discr, vis_order, element_shrink_factor=None)
class meshmode.discretization.visualization.Visualizer(connection, element_shrink_factor=None)
show_scalar_in_mayavi(field, **kwargs)
write_vtk_file(file_name, names_and_fields, compressor=None, real_only=False, overwrite=False)
meshmode.discretization.visualization.write_nodal_adjacency_vtk_file(file_name, mesh, compressor=None, overwrite=False)