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)

Element Groups for 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)

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)
show_scalar_in_matplotlib_3d(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)