Discretizations

Abstract interface

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

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

mesh_el_group
order
index
nelements
nunit_dofs

The number of (for now: nodal) degrees of freedom (“DOFs”) associated with a single element.

ndofs

The total number of (for now: nodal) degrees of freedom (“DOFs”) associated with the element group.

dim
unit_nodes()

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

weights()

Returns an array of length nunit_dofs containing quadrature weights.

is_affine

A bool flag that is True if the local-to-global parametrization of all the elements in the group is affine. Based on meshmode.mesh.MeshElementGroup.is_affine.

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

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

mode_ids()

Return an immutable sequence of opaque (hashable) mode identifiers, one per element of the basis(). The meaning of the mode identifiers is defined by the concrete element group.

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.Discretization(actx: meshmode.array_context.ArrayContext, mesh, group_factory, real_dtype=<class 'numpy.float64'>)

An unstructured composite discretization.

real_dtype
complex_dtype
mesh
dim
ambient_dim
ndofs
groups
empty(actx: meshmode.array_context.ArrayContext, dtype=None)

Return an empty DOFArray.

Parameters

dtype – type special value ‘c’ will result in a vector of dtype complex_dtype. If None (the default), a real vector will be returned.

zeros(actx: meshmode.array_context.ArrayContext, dtype=None)

Return a zero-initialized DOFArray.

Parameters

dtype – type special value ‘c’ will result in a vector of dtype complex_dtype. If None (the default), a real vector will be returned.

empty_like(array: meshmode.dof_array.DOFArray)
zeros_like(array: meshmode.dof_array.DOFArray)
nodes()
Returns

object array of shape (ambient_dim,) containing DOFArrays of node coordinates.

num_reference_derivative(ref_axes, vec)
quad_weights()
Returns

A DOFArray with quadrature weights.

Element Groups for Composite Polynomial Discretization

Group types

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

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.

The mode_ids() are a tuple (one entry per dimension) of directional polynomial degrees on the reference element.

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

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.

The mode_ids() are a tuple (one entry per dimension) of directional polynomial degrees on the reference element.

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

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.

Uses modepy.warp_and_blend_nodes().

The mode_ids() are a tuple (one entry per dimension) of directional polynomial degrees on the reference element.

class meshmode.discretization.poly_element.PolynomialRecursiveNodesElementGroup(mesh_el_group, order, family, index)

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. Depending on the family argument, nodes may be present on the boundary of the simplex. See [Isaac20] for details.

Supports a choice of the base family of 1D nodes, see the documentation of the family argument to recursivenodes.recursive_nodes().

Requires recursivenodes to be installed.

The mode_ids() are a tuple (one entry per dimension) of directional polynomial degrees on the reference element.

Isaac20

Tobin Isaac. Recursive, parameter-free, explicitly defined interpolation nodes for simplices. Arxiv preprint.

New in version 2020.2.

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

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.

The mode_ids() are a tuple (one entry per dimension) of directional polynomial degrees on the reference element.

New in version 2016.1.

class meshmode.discretization.poly_element.LegendreGaussLobattoTensorProductElementGroup(mesh_el_group, order, index)
class meshmode.discretization.poly_element.EquidistantTensorProductElementGroup(mesh_el_group, order, index)

Group factories

class meshmode.discretization.poly_element.ElementGroupFactory
__call__(mesh_ele_group, node_nr_base)
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.PolynomialRecursiveNodesGroupFactory(order, family)
class meshmode.discretization.poly_element.PolynomialEquidistantSimplexGroupFactory(order)

New in version 2016.1.

class meshmode.discretization.poly_element.LegendreGaussLobattoTensorProductGroupFactory(order)

Visualization

meshmode.discretization.visualization.make_visualizer(actx, discr, vis_order, element_shrink_factor=None, force_equidistant=False)
Parameters
  • vis_order – order of the visualization DOFs.

  • element_shrink_factor – number in \((0, 1]\).

  • force_equidistant – if True, the visualization is done on equidistant nodes. If plotting high-order Lagrange VTK elements, this needs to be set to True.

class meshmode.discretization.visualization.Visualizer(connection, element_shrink_factor=None, is_equidistant=False)
show_scalar_in_mayavi(field, **kwargs)
show_scalar_in_matplotlib_3d(field, **kwargs)
write_vtk_file(file_name, names_and_fields, compressor=None, real_only=False, overwrite=False, use_high_order=None, par_manifest_filename=None, par_file_names=None)

Write a Vtk XML file (typical extension .vtu) containing the visualization data in names_and_fields. Can optionally also write manifests for distributed memory simulation (typical extension .pvtu). See also write_parallel_vtk_file() for a convenience wrapper.

Parameters
  • names_and_fields – A list of tuples (name, value), where name is a string and value is a DOFArray or a constant, or an object array of those. value may also be a data class (see dataclasses), whose attributes will be inserted into the visualization with their names prefixed by name. If value is None, then there is no data to write and the corresponding name will not appear in the data file. If value is None, it should be None collectively across all ranks for parallel writes; otherwise the behavior of this routine is undefined.

  • overwrite – If True, silently overwrite existing files.

  • use_high_order – Writes arbitrary order Lagrange VTK elements. These elements are described in this blog post and are available in VTK 8.1 and newer.

  • par_manifest_filename – If not None write a distributed-memory manifest with this file name if file_name matches the first entry in par_file_names.

  • par_file_names – A list of file names of visualization files to include in the distributed-memory manifest.

    Changed in version 2020.2:
  • Added par_manifest_filename and par_file_names.

  • Added use_high_order.

write_parallel_vtk_file(mpi_comm, file_name_pattern, names_and_fields, compressor=None, real_only=False, overwrite=False, use_high_order=None, par_manifest_filename=None)

A convenience wrapper around write_vtk_file() for distributed-memory visualization.

Parameters
  • mpi_comm – An object that supports mpi_comm.Get_rank() and mpi_comm.Get_size() method calls, typically (but not necessarily) an instance of mpi4py.Comm. This is used to determine the current rank as well as the total number of files being written. May also be None in which case a unit-size communicator is assumed.

  • file_name_pattern – A file name pattern (required to end in .vtu) that will be used with str.format() with an (integer) argument of rank to obtain the per-rank file name.

  • par_manifest_filename – as in write_vtk_file(). If not given, par_manifest_filename is synthesized by substituting rank 0 into file_name_pattern and replacing the file extension with .pvtu.

See write_vtk_file() for the meaning of the remainder of the arguments.

New in version 2020.2.

meshmode.discretization.visualization.write_nodal_adjacency_vtk_file(file_name, mesh, compressor=None, overwrite=False)