Common infrastructure

exception meshmode.Error
exception meshmode.DataUnavailable
exception meshmode.FileExistsError

Mesh management

class meshmode.mesh.MeshElementGroup(order, vertex_indices, nodes, element_nr_base=None, node_nr_base=None, unit_nodes=None, dim=None)

A group of elements sharing a common reference element.

order
vertex_indices

An array (nelements, ref_element.nvertices) of (mesh-wide) vertex indices.

nodes

An array of node coordinates with shape (mesh.ambient_dim, nelements, nunit_nodes).

unit_nodes

(dim, nunit_nodes)

element_nr_base

Lowest element number in this element group.

node_nr_base

Lowest node number in this element group.

dim

The number of dimensions spanned by the element. Not the ambient dimension, see Mesh.ambient_dim for that.

face_vertex_indices()

Return a tuple of tuples indicating which vertices (in mathematically positive ordering) make up each face of an element in this group.

vertex_unit_coordinates()

Return an array of shape (nfaces, dim) with the unit coordinates of each vertex.

nfaces
__eq__(other)

Return self==value.

__ne__(other)

Return self!=value.

copy(**kwargs)
dim
face_vertex_indices()

Return a tuple of tuples indicating which vertices (in mathematically positive ordering) make up each face of an element in this group.

join_mesh(element_nr_base, node_nr_base)
nelements
nfaces
nnodes
nunit_nodes
vertex_unit_coordinates()

Return an array of shape (nfaces, dim) with the unit coordinates of each vertex.

class meshmode.mesh.Mesh(vertices, groups, skip_tests=False, node_vertex_consistency_tolerance=None, skip_element_orientation_test=False, nodal_adjacency=None, facial_adjacency_groups=None, boundary_tags=None, vertex_id_dtype=<class 'numpy.int32'>, element_id_dtype=<class 'numpy.int32'>, is_conforming=None)
vertices

An array of vertex coordinates with shape (ambient_dim, nvertices)

groups

A list of MeshElementGroup instances.

nodal_adjacency

An instance of NodalAdjacency.

Referencing this attribute may raise meshmode.DataUnavailable.

facial_adjacency_groups

A list of mappings from neighbor group IDs to instances of FacialAdjacencyGroup.

facial_adjacency_groups[igrp][ineighbor_group] gives the set of facial adjacency relations between group igrp and ineighbor_group. ineighbor_group and igrp may be identical, or ineighbor_group may be None, in which case an :class:InterPartitionAdjacency group containing boundary faces is returned.

Referencing this attribute may raise meshmode.DataUnavailable.

_images/facial-adjacency-group.png

For example for the mesh in the figure, the following data structure would be present:

[
    {...},  # connectivity for group 0
    {...},  # connectivity for group 1
    {...},  # connectivity for group 2
    {       # connectivity for group 3
        1: FacialAdjacencyGroup(...)  # towards group 1, green
        2: FacialAdjacencyGroup(...)  # towards group 2, pink
        None: FacialAdjacencyGroup(...)  # towards the boundary, orange
    }
]

(Note that element groups are not necessarily contiguous like the figure may suggest.)

boundary_tags

A tuple of boundary tag identifiers. BTAG_ALL and BTAG_REALLY_ALL are guranateed to exist.

btag_to_index

A mapping that maps boundary tag identifiers to their corresponding index.

vertex_id_dtype
element_id_dtype
is_conforming

True if it is known that all element interfaces are conforming. False if it is known that some element interfaces are non-conforming. None if neither of the two is known.

__eq__(other)

Return self==value.

__ne__(other)

Return self!=value.

ambient_dim
boundary_tag_bit(boundary_tag)
dim
face_id_dtype

alias of numpy.int8

facial_adjacency_groups
get_copy_kwargs(**kwargs)
nelements
nodal_adjacency
nodal_adjacency_init_arg()

Returns an ‘nodal_adjacency’ argument that can be passed to a Mesh constructor.

nvertices
class meshmode.mesh.NodalAdjacency(valuedict=None, exclude=['self'], **kwargs)

Describes nodal element adjacency information, i.e. information about elements that touch in at least one point.

neighbors_starts

element_id_t [nelements+1]

Use together with neighbors. neighbors_starts[iel] and neighbors_starts[iel+1] together indicate a ranges of element indices neighbors which are adjacent to iel.

neighbors

element_id_t []

See neighbors_starts.

__eq__(other)

Return self==value.

__ne__(other)

Return self!=value.

class meshmode.mesh.FacialAdjacencyGroup(valuedict=None, exclude=['self'], **kwargs)

Describes facial element adjacency information for one MeshElementGroup, i.e. information about elements that share (part of) a face.

_images/facial-adjacency-group.png

Represents (for example) one of the (colored) interfaces between MeshElementGroup instances, or an interface between MeshElementGroup and a boundary. (Note that element groups are not necessarily contiguous like the figure may suggest.)

igroup
ineighbor_group

ID of neighboring group, or None for boundary faces. If identical to igroup, then this contains the self-connectivity in this group.

elements

element_id_t [nfagrp_elements]. elements[i] Group-local element numbers.

element_faces

face_id_t [nfagrp_elements]. element_faces[i] indicate what face index of the opposite element indicated in neighbors[iface] touches face number iface of element number iel_grp in this element group.

neighbors

element_id_t [nfagrp_elements]. neighbors[i] gives the element number within ineighbor_group of the element opposite elements[i].

If this number is negative, then this indicates that this is a boundary face, and the bits set in -neighbors[i] should be interpreted according to Mesh.boundary_tags.

neighbor_faces

face_id_t [nfagrp_elements]. neighbor_faces[i] indicate what face index of the opposite element indicated in neighbors[i] has facial contact with face number element_faces[i] of element number elements[i] in element group igroup.

Zero if neighbors[i] is negative.

__eq__(other)

Return self==value.

__ne__(other)

Return self!=value.

meshmode.mesh.as_python(mesh, function_name='make_mesh')

Return a snippet of Python code (as a string) that will recreate the mesh given as an input parameter.

meshmode.mesh.check_bc_coverage(mesh, boundary_tags, incomplete_ok=False, true_boundary_only=True)

Verify boundary condition coverage.

Given a list of boundary tags as boundary_tags, this function verifies that

1. the union of all these boundaries gives the complete boundary, 1. all these boundaries are disjoint.
Parameters:
  • incomplete_ok – Do not report an error if some faces are not covered by the boundary conditions.
  • true_boundary_only – only verify for faces tagged with BTAG_ALL.
meshmode.mesh.is_boundary_tag_empty(mesh, boundary_tag)

Return True if the corresponding boundary tag does not occur as part of mesh.

Predefined Boundary tags

class meshmode.mesh.BTAG_NONE

A boundary tag representing an empty boundary or volume.

class meshmode.mesh.BTAG_ALL

A boundary tag representing the entire boundary or volume.

In the case of the boundary, TAG_ALL does not include rank boundaries, or, more generally, anything tagged with TAG_NO_BOUNDARY.

class meshmode.mesh.BTAG_REALLY_ALL

A boundary tag representing the entire boundary.

Unlike TAG_ALL, this includes rank boundaries, or, more generally, everything tagged with TAG_NO_BOUNDARY.

class meshmode.mesh.BTAG_NO_BOUNDARY

A boundary tag indicating that this edge should not fall under TAG_ALL. Among other things, this is used to keep rank boundaries out of BTAG_ALL.

class meshmode.mesh.BTAG_PARTITION(part_nr)

A boundary tag indicating that this edge is adjacent to an element of another Mesh. The partition number of the adjacent mesh is given by part_nr.

part_nr

New in version 2017.1.

Design of the Data Structure

Why does a Mesh need to be broken into MeshElementGroup instances?

Elements can be of different types (e.g. triangle, quadrilateral, tetrahedron, what have you). In addition, elements may vary in the polynomial degree used to represent them (see also below).

All these bits of information could in principle be stored by element, but having large, internally homogeneous groups is a good thing from an efficiency standpoint. (So that you can, e.g., launch one GPU kernel to deal with all order-3 triangles, instead of maybe having to dispatch based on type and size inside the kernel)

What is the difference between ‘vertices’ and ‘nodes’?

Nodes exist mainly to represent the (potentially non-affine) deformation of each element, by a one-to-one correspondence with MeshElementGroup.unit_nodes. They are unique to each element. Vertices on the other hand exist to clarify whether or not a point shared by two elements is actually identical (or just happens to be “close”). This is done by assigning (single, globally shared) vertex numbers and having elements refer to them.

Consider the following picture:

_images/nodes-vertices.png

Mesh Data Structure

class meshmode.mesh.MeshElementGroup(order, vertex_indices, nodes, element_nr_base=None, node_nr_base=None, unit_nodes=None, dim=None)

A group of elements sharing a common reference element.

order
vertex_indices

An array (nelements, ref_element.nvertices) of (mesh-wide) vertex indices.

nodes

An array of node coordinates with shape (mesh.ambient_dim, nelements, nunit_nodes).

unit_nodes

(dim, nunit_nodes)

element_nr_base

Lowest element number in this element group.

node_nr_base

Lowest node number in this element group.

dim

The number of dimensions spanned by the element. Not the ambient dimension, see Mesh.ambient_dim for that.

face_vertex_indices()

Return a tuple of tuples indicating which vertices (in mathematically positive ordering) make up each face of an element in this group.

vertex_unit_coordinates()

Return an array of shape (nfaces, dim) with the unit coordinates of each vertex.

nfaces
__eq__(other)

Return self==value.

__ne__(other)

Return self!=value.

copy(**kwargs)
dim
face_vertex_indices()

Return a tuple of tuples indicating which vertices (in mathematically positive ordering) make up each face of an element in this group.

join_mesh(element_nr_base, node_nr_base)
nelements
nfaces
nnodes
nunit_nodes
vertex_unit_coordinates()

Return an array of shape (nfaces, dim) with the unit coordinates of each vertex.

class meshmode.mesh.Mesh(vertices, groups, skip_tests=False, node_vertex_consistency_tolerance=None, skip_element_orientation_test=False, nodal_adjacency=None, facial_adjacency_groups=None, boundary_tags=None, vertex_id_dtype=<class 'numpy.int32'>, element_id_dtype=<class 'numpy.int32'>, is_conforming=None)
vertices

An array of vertex coordinates with shape (ambient_dim, nvertices)

groups

A list of MeshElementGroup instances.

nodal_adjacency

An instance of NodalAdjacency.

Referencing this attribute may raise meshmode.DataUnavailable.

facial_adjacency_groups

A list of mappings from neighbor group IDs to instances of FacialAdjacencyGroup.

facial_adjacency_groups[igrp][ineighbor_group] gives the set of facial adjacency relations between group igrp and ineighbor_group. ineighbor_group and igrp may be identical, or ineighbor_group may be None, in which case an :class:InterPartitionAdjacency group containing boundary faces is returned.

Referencing this attribute may raise meshmode.DataUnavailable.

_images/facial-adjacency-group.png

For example for the mesh in the figure, the following data structure would be present:

[
    {...},  # connectivity for group 0
    {...},  # connectivity for group 1
    {...},  # connectivity for group 2
    {       # connectivity for group 3
        1: FacialAdjacencyGroup(...)  # towards group 1, green
        2: FacialAdjacencyGroup(...)  # towards group 2, pink
        None: FacialAdjacencyGroup(...)  # towards the boundary, orange
    }
]

(Note that element groups are not necessarily contiguous like the figure may suggest.)

boundary_tags

A tuple of boundary tag identifiers. BTAG_ALL and BTAG_REALLY_ALL are guranateed to exist.

btag_to_index

A mapping that maps boundary tag identifiers to their corresponding index.

vertex_id_dtype
element_id_dtype
is_conforming

True if it is known that all element interfaces are conforming. False if it is known that some element interfaces are non-conforming. None if neither of the two is known.

__eq__(other)

Return self==value.

__ne__(other)

Return self!=value.

ambient_dim
boundary_tag_bit(boundary_tag)
dim
face_id_dtype

alias of numpy.int8

facial_adjacency_groups
get_copy_kwargs(**kwargs)
nelements
nodal_adjacency
nodal_adjacency_init_arg()

Returns an ‘nodal_adjacency’ argument that can be passed to a Mesh constructor.

nvertices

Mesh generation

Curves

meshmode.mesh.generation.make_curve_mesh(curve_f, element_boundaries, order, unit_nodes=None, node_vertex_consistency_tolerance=None, return_parametrization_points=False)
Parameters:
  • curve_f – A callable representing a parametrization for a curve, accepting a vector of point locations and returning an array of shape (2, npoints).
  • element_boundaries – a vector of element boundary locations in \([0,1]\), in order. 0 must be the first entry, 1 the last one.
  • unit_nodes – if given, the unit nodes to use. Must have shape (dim, nnoodes).
Returns:

a meshmode.mesh.Mesh, or if return_parametrization_points is True, a tuple (mesh, par_points), where par_points is an array of parametrization points.

Curve parametrizations

meshmode.mesh.generation.ellipse(aspect_ratio, t)
Parameters:t – the parametrization, runs from [0,1)
Returns:an array of shape (2, npoints)
meshmode.mesh.generation.cloverleaf(t)
Parameters:t – the parametrization, runs from [0,1)
Returns:an array of shape (2, npoints)
meshmode.mesh.generation.starfish
meshmode.mesh.generation.drop(t)
Parameters:t – the parametrization, runs from [0,1)
Returns:an array of shape (2, npoints)
meshmode.mesh.generation.n_gon(n_corners, t)
Parameters:t – the parametrization, runs from [0,1)
Returns:an array of shape (2, npoints)
meshmode.mesh.generation.qbx_peanut(t)
meshmode.mesh.generation.apple(a, t)
Parameters:
  • a – 0 <= a <= 1/2; roundedness: 0 returns a circle, 1/2 returns a cardioid
  • t – the parametrization, runs from [0,1)
Returns:

an array of shape (2, npoints)

class meshmode.mesh.generation.WobblyCircle(coeffs)
static random(ncoeffs, seed)
__call__(t)
Parameters:t – the parametrization, runs from [0,1)
Returns:an array of shape (2, npoints)

Surfaces

meshmode.mesh.generation.generate_icosahedron(r, order)
meshmode.mesh.generation.generate_icosphere(r, order)
meshmode.mesh.generation.generate_torus(r_outer, r_inner, n_outer=20, n_inner=10, order=1)
meshmode.mesh.generation.refine_mesh_and_get_urchin_warper(order, m, n, est_rel_interp_tolerance, min_rad=0.2, uniform_refinement_rounds=0)
Returns:

a tuple (refiner, warp_mesh), where refiner is a meshmode.refinement.Refiner (from which the unwarped mesh may be obtained), and whose meshmode.refinement.Refiner.get_current_mesh() returns a locally-refined meshmode.mesh.Mesh of a sphere and warp_mesh is a callable taking and returning a mesh that warps the unwarped mesh into a smooth shape govered by a spherical harmonic of order (m, n).

Parameters:
  • order – the polynomial order of the returned mesh
  • est_rel_interp_tolerance – a tolerance for the relative interpolation error estimates on the warped version of the mesh.
meshmode.mesh.generation.generate_urchin(order, m, n, est_rel_interp_tolerance, min_rad=0.2)
Returns:

a refined meshmode.mesh.Mesh of a smooth shape govered by a spherical harmonic of order (m, n).

Parameters:
  • order – the polynomial order of the returned mesh
  • est_rel_interp_tolerance – a tolerance for the relative interpolation error estimates on the warped version of the mesh.

Volumes

meshmode.mesh.generation.generate_box_mesh(axis_coords, order=1, coord_dtype=<class 'numpy.float64'>, group_factory=None)

Create a semi-structured mesh.

Parameters:
  • axis_coords – a tuple with a number of entries corresponding to the number of dimensions, with each entry a numpy array specifying the coordinates to be used along that axis.
  • group_factory – One of meshmode.mesh.SimplexElementGroup or meshmode.mesh.TensorProductElementGroup.

Changed in version 2017.1: group_factory parameter added.

meshmode.mesh.generation.generate_regular_rect_mesh(a=(0, 0), b=(1, 1), n=(5, 5), order=1)

Create a semi-structured rectangular mesh.

Parameters:
  • a – the lower left hand point of the rectangle
  • b – the upper right hand point of the rectangle
  • n – a tuple of integers indicating the total number of points on [a,b].
meshmode.mesh.generation.generate_warped_rect_mesh(dim, order, n)

Generate a mesh of a warped line/square/cube. Mainly useful for testing functionality with curvilinear meshes.

Tools for Iterative Refinement

meshmode.mesh.generation.warp_and_refine_until_resolved(unwarped_mesh_or_refiner, warp_callable, est_rel_interp_tolerance)

Given an original (“un-warped”) meshmode.mesh.Mesh and a warping function warp_callable that takes and returns a mesh and a tolerance to which the mesh should be resolved by the mapping polynomials, this function will iteratively refine the unwarped_mesh until relative interpolation error estimates on the warped version are smaller than est_rel_interp_tolerance on each element.

Returns:The refined, un-warped mesh.

New in version 2018.1.

Mesh input/output

class meshmode.mesh.io.ScriptSource(source, extension)

New in version 2016.1.

class meshmode.mesh.io.FileSource(filename)

New in version 2014.1.

class meshmode.mesh.io.ScriptWithFilesSource(source, filenames, source_name='temp.geo')

New in version 2016.1.

source

The script code to be fed to gmsh.

filenames

The names of files to be copied to the temporary directory where gmsh is run.

meshmode.mesh.io.read_gmsh(filename, force_ambient_dim=None, mesh_construction_kwargs=None)

Read a gmsh mesh file from filename and return a meshmode.mesh.Mesh.

Parameters:
  • force_ambient_dim – if not None, truncate point coordinates to this many dimensions.
  • mesh_construction_kwargsNone or a dictionary of keyword arguments passed to the meshmode.mesh.Mesh constructor.
meshmode.mesh.io.generate_gmsh(source, dimensions=None, order=None, other_options=[], extension='geo', gmsh_executable='gmsh', force_ambient_dim=None, output_file_name='output.msh', mesh_construction_kwargs=None)

Run gmsh on the input given by source, and return a meshmode.mesh.Mesh based on the result.

Parameters:
  • source – an instance of either FileSource or LiteralSource
  • force_ambient_dim – if not None, truncate point coordinates to this many dimensions.
  • mesh_construction_kwargsNone or a dictionary of keyword arguments passed to the meshmode.mesh.Mesh constructor.
meshmode.mesh.io.from_meshpy(mesh_info, order=1)

Imports a mesh from a meshpy mesh_info data structure, which may be generated by either meshpy.triangle or meshpy_tet.

meshmode.mesh.io.from_vertices_and_simplices(vertices, simplices, order=1, fix_orientation=False)

Imports a mesh from a numpy array of vertices and an array of simplices.

Parameters:
  • vertices – An array of vertex coordinates with shape (ambient_dim, nvertices)
  • simplices – An array (nelements, nvertices) of (mesh-wide) vertex indices.
meshmode.mesh.io.to_json(mesh)

Return a JSON-able Python data structure for mesh. The structure directly reflects the Mesh data structure.

Mesh processing

meshmode.mesh.processing.find_group_indices(groups, meshwide_elems)
Parameters:
  • groups – A list of :class:MeshElementGroup instances that contain meshwide_elems.
  • meshwide_elems – A :class:numpy.ndarray of mesh-wide element numbers Usually computed by elem + element_nr_base.
Returns:

A :class:numpy.ndarray of group numbers that meshwide_elem belongs to.

meshmode.mesh.processing.partition_mesh(mesh, part_per_element, part_num)
Parameters:
  • mesh – A meshmode.mesh.Mesh to be partitioned.
  • part_per_element – A numpy.ndarray containing one integer per element of mesh indicating which part of the partitioned mesh the element is to become a part of.
  • part_num – The part number of the mesh to return.
Returns:

A tuple (part_mesh, part_to_global), where part_mesh is a meshmode.mesh.Mesh that is a partition of mesh, and part_to_global is a numpy.ndarray mapping element numbers on part_mesh to ones in mesh.

New in version 2017.1.

meshmode.mesh.processing.find_volume_mesh_element_orientations(mesh, tolerate_unimplemented_checks=False)

Return a positive floating point number for each positively oriented element, and a negative floating point number for each negatively oriented element.

Parameters:tolerate_unimplemented_checks – If True, elements for which no check is available will return NaN.
meshmode.mesh.processing.perform_flips(mesh, flip_flags, skip_tests=False)
Parameters:flip_flags – A numpy.ndarray with mesh.nelements entries indicating by their Boolean value whether the element is to be flipped.
meshmode.mesh.processing.find_bounding_box(mesh)
Returns:a tuple (min, max), each consisting of a numpy.ndarray indicating the minimal and maximal extent of the geometry along each axis.
meshmode.mesh.processing.merge_disjoint_meshes(meshes, skip_tests=False, single_group=False)
meshmode.mesh.processing.map_mesh(mesh, f)

Apply the map f to the mesh. f needs to accept and return arrays of shape (ambient_dim, npoints).

meshmode.mesh.processing.affine_map(mesh, A=None, b=None)

Apply the affine map f(x)=Ax+b to the geometry of mesh.

Mesh refinement

class meshmode.mesh.refinement.Refiner(mesh)
class meshmode.mesh.refinement.RefinerWithoutAdjacency(mesh)

A refiner that may be applied to non-conforming meshmode.mesh.Mesh instances. It does not generate adjacency information, and it is typically faster than meshmode.mesh.refinement.Refiner.

Note

If the input meshes to this refiner are not conforming, then the resulting meshes may contain duplicated vertices. (I.e. two different numbers referring to the same geometric vertex.)

meshmode.mesh.refinement.refine_uniformly(mesh, iterations)

Mesh visualization

meshmode.mesh.visualization.draw_2d_mesh(mesh, draw_vertex_numbers=True, draw_element_numbers=True, draw_nodal_adjacency=False, draw_face_numbers=False, set_bounding_box=False, **kwargs)
meshmode.mesh.visualization.draw_curve(mesh, el_bdry_style='o', el_bdry_kwargs=None, node_style='x-', node_kwargs=None)
meshmode.mesh.visualization.write_vertex_vtk_file(mesh, file_name, compressor=None, overwrite=False)
meshmode.mesh.visualization.mesh_to_tikz(mesh)