Mesh Classes¶
This page documents the mesh classes provided by mmgpy.
Unified Mesh Class¶
The Mesh class is the primary public API for mmgpy. It provides a unified interface that auto-detects mesh type based on input data:
mmgpy.Mesh
¶
Mesh(
source: NDArray[floating] | str | Path | UnstructuredGrid | PolyData,
cells: NDArray[integer] | None = None,
)
Unified mesh class with auto-detection of mesh type.
This class provides a single interface for working with 2D planar, 3D volumetric, and 3D surface meshes. The mesh type is automatically detected from the input data.
Parameters¶
source : ndarray | str | Path | pv.UnstructuredGrid | pv.PolyData
Either:
- Vertex coordinates array (requires cells parameter)
- File path to load mesh from
- PyVista mesh object
cells : ndarray, optional
Cell connectivity array. Required when source is vertices.
Attributes¶
kind : MeshKind The type of mesh (TETRAHEDRAL, TRIANGULAR_2D, or TRIANGULAR_SURFACE).
Examples¶
Create a mesh from vertices and cells:
vertices = np.array([[0, 0, 0], [1, 0, 0], [0, 1, 0], [0, 0, 1]]) cells = np.array([[0, 1, 2, 3]]) mesh = Mesh(vertices, cells) mesh.kind MeshKind.TETRAHEDRAL
Load a mesh from file:
mesh = Mesh("mesh.vtk")
Create from PyVista:
pv_mesh = pv.read("mesh.vtk") mesh = Mesh(pv_mesh)
Initialize a Mesh from various sources.
vtk
property
¶
Get the PyVista mesh representation.
This property provides direct access to the PyVista mesh for use with custom plotters or other PyVista operations.
Returns¶
pv.UnstructuredGrid | pv.PolyData PyVista mesh object.
Examples¶
plotter = pv.Plotter() plotter.add_mesh(mesh.vtk, show_edges=True) plotter.show()
get_vertices
¶
get_vertices_with_refs
¶
Get vertex coordinates and reference markers.
Returns¶
vertices : ndarray Vertex coordinates. refs : ndarray Reference markers for each vertex.
set_vertices
¶
Set vertex coordinates.
Parameters¶
vertices : ndarray Vertex coordinates. refs : ndarray, optional Reference markers for each vertex.
get_triangles
¶
get_triangles_with_refs
¶
Get triangle connectivity and reference markers.
Returns¶
triangles : ndarray Triangle connectivity. refs : ndarray Reference markers for each triangle.
set_triangles
¶
Set triangle connectivity.
Parameters¶
triangles : ndarray Triangle connectivity (Nx3). refs : ndarray, optional Reference markers for each triangle.
get_edges
¶
get_edges_with_refs
¶
Get edge connectivity and reference markers.
Returns¶
edges : ndarray Edge connectivity. refs : ndarray Reference markers for each edge.
set_edges
¶
Set edge connectivity.
Parameters¶
edges : ndarray Edge connectivity (Nx2). refs : ndarray, optional Reference markers for each edge.
get_tetrahedra
¶
get_tetrahedra_with_refs
¶
get_elements
¶
get_elements_with_refs
¶
set_field
¶
Set a solution field.
Parameters¶
key : str Field name. value : ndarray Field values (one per vertex).
get_field
¶
__setitem__
¶
Set a solution field using dictionary syntax.
__getitem__
¶
Get a solution field using dictionary syntax.
set_user_field
¶
Set a user-defined field for transfer during remeshing.
Unlike MMG's built-in fields (metric, displacement, levelset), user fields are arbitrary data arrays that can be transferred to the new mesh after remeshing via interpolation.
Parameters¶
name : str Field name (any string except reserved names like "metric"). values : ndarray Field values, shape (n_vertices,) for scalars or (n_vertices, n_components) for vectors/tensors.
Examples¶
mesh.set_user_field("temperature", temperature_array) mesh.set_user_field("velocity", velocity_array) # (N, 3) vector mesh.remesh(hmax=0.1, transfer_fields=True) new_temp = mesh.get_user_field("temperature")
get_user_field
¶
get_user_fields
¶
has_user_field
¶
get_bounds
¶
Get the bounding box of the mesh.
Returns¶
tuple[ndarray, ndarray] Tuple of (min_coords, max_coords), each shape (3,) for 3D or (2,) for 2D meshes.
Raises¶
ValueError If the mesh has no vertices.
Examples¶
mesh = Mesh(vertices, cells) min_pt, max_pt = mesh.get_bounds() print(f"Size: {max_pt - min_pt}")
get_center_of_mass
¶
Get the centroid (center of mass) of the mesh.
For volume meshes, computes volume-weighted centroid. For surface meshes, computes area-weighted centroid. For 2D meshes, computes area-weighted centroid.
Returns¶
ndarray Centroid coordinates, shape (3,) for 3D or (2,) for 2D.
Examples¶
mesh = Mesh(vertices, cells) center = mesh.get_center_of_mass() print(f"Center: {center}")
compute_volume
¶
compute_surface_area
¶
Compute the total surface area of the mesh.
For volume meshes (TETRAHEDRAL), computes boundary surface area (triangles stored in the mesh represent boundary faces). For surface meshes (TRIANGULAR_SURFACE), computes total area. For 2D meshes (TRIANGULAR_2D), computes total area.
Returns¶
float Total surface area in mesh units squared.
Examples¶
mesh = Mesh(vertices, cells) area = mesh.compute_surface_area() print(f"Surface area: {area:.2f} mm^2")
get_diagonal
¶
get_adjacent_elements
¶
get_vertex_neighbors
¶
get_element_quality
¶
get_element_qualities
¶
save
¶
Save mesh to file.
Parameters¶
filename : str or Path Output file path. Format determined by extension.
remesh
¶
remesh(
options: Mmg3DOptions | Mmg2DOptions | MmgSOptions | None = None,
*,
progress: ProgressParam = True,
transfer_fields: FieldTransferParam = False,
interpolation: str = "linear",
**kwargs: Any,
) -> RemeshResult
Remesh the mesh in-place.
Parameters¶
options : Mmg3DOptions | Mmg2DOptions | MmgSOptions, optional Options object for remeshing parameters. progress : bool | Callable[[ProgressEvent], bool] | None, default=True Progress reporting option: - True: Show Rich progress bar (default) - False or None: No progress reporting - Callable: Custom callback that receives ProgressEvent and returns True to continue or False to cancel transfer_fields : bool | Sequence[str] | None, default=False Transfer user-defined fields to the new mesh via interpolation: - False or None: No field transfer (default), clears existing user fields - True: Transfer all user fields - List of field names: Transfer only specified fields interpolation : str, default="linear" Interpolation method for field transfer: - "linear": Barycentric interpolation (recommended) - "nearest": Nearest vertex value **kwargs : float Individual remeshing parameters (hmin, hmax, hsiz, hausd, etc.).
Notes¶
Memory: When transfer_fields is enabled, the original mesh vertices
and elements are copied before remeshing, temporarily doubling memory
usage for large meshes.
Surface meshes (TRIANGULAR_SURFACE): Field transfer uses 3D Delaunay triangulation for point location, which may not work well for nearly planar surface meshes. Consider using volumetric meshes for field transfer.
Returns¶
RemeshResult Statistics from the remeshing operation.
Raises¶
CancellationError If the progress callback returns False to cancel the operation.
Examples¶
mesh.remesh(hmax=0.1) # Shows progress bar by default
mesh.remesh(hmax=0.1, progress=False) # No progress bar
Transfer fields during remeshing¶
mesh.set_user_field("temperature", temperature_array) mesh.remesh(hmax=0.1, transfer_fields=True) new_temp = mesh.get_user_field("temperature")
Transfer specific fields¶
mesh.remesh(hmax=0.1, transfer_fields=["temperature", "velocity"])
def my_callback(event): ... print(f"{event.phase}: {event.message}") ... return True # Continue mesh.remesh(hmax=0.1, progress=my_callback)
remesh_lagrangian
¶
remesh_lagrangian(
displacement: NDArray[float64], *, progress: ProgressParam = True, **kwargs: Any
) -> RemeshResult
Remesh with Lagrangian motion.
Only available for TETRAHEDRAL and TRIANGULAR_2D meshes.
Parameters¶
displacement : ndarray Displacement field for each vertex. progress : bool | Callable[[ProgressEvent], bool] | None, default=True Progress reporting option: - True: Show Rich progress bar (default) - False or None: No progress reporting - Callable: Custom callback that receives ProgressEvent and returns True to continue or False to cancel **kwargs : float Additional remeshing parameters.
Returns¶
RemeshResult Statistics from the remeshing operation.
Raises¶
TypeError If mesh is TRIANGULAR_SURFACE. CancellationError If the progress callback returns False to cancel the operation.
remesh_levelset
¶
remesh_levelset(
levelset: NDArray[float64], *, progress: ProgressParam = True, **kwargs: Any
) -> RemeshResult
Remesh with level-set discretization.
Parameters¶
levelset : ndarray Level-set field for each vertex. progress : bool | Callable[[ProgressEvent], bool] | None, default=True Progress reporting option: - True: Show Rich progress bar (default) - False or None: No progress reporting - Callable: Custom callback that receives ProgressEvent and returns True to continue or False to cancel **kwargs : float Additional remeshing parameters.
Returns¶
RemeshResult Statistics from the remeshing operation.
Raises¶
CancellationError If the progress callback returns False to cancel the operation.
remesh_optimize
¶
Optimize mesh quality without changing topology.
Only moves vertices to improve element quality. No points are inserted or removed.
Parameters¶
progress : bool | Callable[[ProgressEvent], bool] | None, default=True Progress reporting option: - True: Show Rich progress bar (default) - False or None: No progress reporting - Callable: Custom callback verbose : int | None Verbosity level (-1=silent, 0=errors, 1=info).
Returns¶
RemeshResult Statistics from the remeshing operation.
remesh_uniform
¶
remesh_uniform(
size: float, *, progress: ProgressParam = True, verbose: int | None = None
) -> RemeshResult
Remesh with uniform element size.
Parameters¶
size : float Target edge size for all elements. progress : bool | Callable[[ProgressEvent], bool] | None, default=True Progress reporting option: - True: Show Rich progress bar (default) - False or None: No progress reporting - Callable: Custom callback verbose : int | None Verbosity level (-1=silent, 0=errors, 1=info).
Returns¶
RemeshResult Statistics from the remeshing operation.
set_size_sphere
¶
Set target edge size within a spherical region.
Parameters¶
center : array-like Center point of the sphere. radius : float Radius of the sphere. size : float Target edge size within the sphere.
set_size_box
¶
Set target edge size within a box region.
Parameters¶
bounds : array-like Bounding box as [[xmin, ymin, zmin], [xmax, ymax, zmax]]. size : float Target edge size within the box.
set_size_cylinder
¶
set_size_cylinder(
point1: Sequence[float] | NDArray[float64],
point2: Sequence[float] | NDArray[float64],
radius: float,
size: float,
) -> None
Set target edge size within a cylindrical region.
Only available for TETRAHEDRAL and TRIANGULAR_SURFACE meshes.
Parameters¶
point1 : array-like First endpoint of the cylinder axis. point2 : array-like Second endpoint of the cylinder axis. radius : float Radius of the cylinder. size : float Target edge size within the cylinder.
Raises¶
TypeError If mesh is TRIANGULAR_2D.
set_size_from_point
¶
set_size_from_point(
point: Sequence[float] | NDArray[float64],
near_size: float,
far_size: float,
influence_radius: float,
) -> None
Set target edge size based on distance from a point.
Parameters¶
point : array-like Reference point. near_size : float Target edge size at the reference point. far_size : float Target edge size at the influence radius. influence_radius : float Radius of influence.
get_local_sizing_count
¶
apply_local_sizing
¶
Apply local sizing constraints to the metric field.
This is called automatically before remeshing if sizing constraints have been added.
to_pyvista
¶
plot
¶
validate
¶
validate(
*,
detailed: bool = False,
strict: bool = False,
check_geometry: bool = True,
check_topology: bool = True,
check_quality: bool = True,
min_quality: float = 0.1,
) -> bool | ValidationReport
Validate the mesh and check for issues.
Parameters¶
detailed : bool If True, return a ValidationReport with detailed information. If False, return a simple boolean. strict : bool If True, raise ValidationError on any issue (including warnings). check_geometry : bool Check for geometric issues (inverted/degenerate elements). check_topology : bool Check for topological issues (orphan vertices, non-manifold edges). check_quality : bool Check element quality against threshold. min_quality : float Minimum acceptable element quality (0-1).
Returns¶
bool | ValidationReport If detailed=False, returns True if valid, False otherwise. If detailed=True, returns full ValidationReport.
Raises¶
ValidationError If strict=True and any issues are found.
Examples¶
mesh = Mesh(vertices, cells) if mesh.validate(): ... print("Mesh is valid")
report = mesh.validate(detailed=True) print(f"Quality: {report.quality.mean:.3f}")
edit_sizing
¶
edit_sizing(
*, mode: str = "sphere", default_size: float = 0.01, default_radius: float = 0.1
) -> None
Launch interactive sizing editor.
Opens a PyVista window for visually defining local sizing constraints by clicking on the mesh.
Parameters¶
mode : str Initial interaction mode: "sphere", "box", "cylinder", or "point". default_size : float Default target edge size for constraints. default_radius : float Default radius for sphere and cylinder constraints.
Examples¶
mesh = Mesh(vertices, cells) mesh.edit_sizing() # Opens interactive editor mesh.remesh() # Uses interactively defined sizing
__enter__
¶
__exit__
¶
__exit__(
exc_type: type[BaseException] | None,
exc_val: BaseException | None,
exc_tb: TracebackType | None,
) -> bool
Exit the context manager.
Currently performs no cleanup, but provides a consistent API for resource management patterns.
Returns¶
bool False, to not suppress any exceptions.
checkpoint
¶
Create a checkpoint for transactional modifications.
Returns a context manager that captures the current mesh state.
On exit, if commit() was not called or an exception occurred,
the mesh is automatically rolled back to its checkpoint state.
Returns¶
MeshCheckpoint A context manager for transactional mesh modifications.
Notes¶
The checkpoint stores a complete copy of the mesh data including vertices, elements, reference markers, and solution fields (metric, displacement, levelset). For large meshes, this may consume significant memory.
Note: The tensor field is not saved because it shares memory with metric in MMG's internal representation.
Examples¶
mesh = Mesh(vertices, cells) with mesh.checkpoint() as snapshot: ... mesh.remesh(hmax=0.01) ... if mesh.validate(): ... snapshot.commit() # Keep changes ... # If not committed, changes are rolled back
Automatic rollback on exception¶
with mesh.checkpoint(): ... mesh.remesh(hmax=0.01) ... raise ValueError("Simulated failure")
mesh is restored to original state¶
copy
¶
Create a working copy that is discarded on exit.
Returns a context manager that yields a copy of the mesh.
The copy can be freely modified without affecting the original.
Use update_from() to apply changes from the copy to the original.
Yields¶
Mesh A copy of this mesh.
Examples¶
original = Mesh(vertices, cells) with original.copy() as working: ... working.remesh(hmax=0.1) ... if len(working.get_vertices()) < len(original.get_vertices()) * 2: ... original.update_from(working)
working is discarded on exit¶
update_from
¶
Update this mesh from another mesh's state.
Replaces the vertices and elements of this mesh with those from the other mesh. Both meshes must be of the same kind.
Parameters¶
other : Mesh The mesh to copy state from.
Raises¶
TypeError If the meshes are of different kinds.
Examples¶
original = Mesh(vertices, cells) with original.copy() as working: ... working.remesh(hmax=0.1) ... original.update_from(working)
options: members: - init - kind - save - get_vertices - get_triangles - get_tetrahedra - get_edges - get_mesh_size - remesh - remesh_optimize - remesh_uniform - remesh_levelset - remesh_lagrangian - validate - to_pyvista - set_size_sphere - set_size_box - set_size_cylinder - set_size_from_point - clear_local_sizing - get_local_sizing_count
Mesh Kind Enumeration¶
mmgpy.MeshKind
¶
Bases: Enum
Enumeration of mesh types.
Attributes¶
TETRAHEDRAL 3D volumetric mesh with tetrahedral elements. TRIANGULAR_2D 2D planar mesh with triangular elements. TRIANGULAR_SURFACE 3D surface mesh with triangular elements.
options: show_root_heading: true
Usage Examples¶
Creating Meshes¶
import mmgpy
import numpy as np
# From file
mesh = mmgpy.Mesh("input.mesh")
# From arrays - 3D tetrahedral mesh
vertices = np.array([
[0, 0, 0],
[1, 0, 0],
[0, 1, 0],
[0, 0, 1],
], dtype=np.float64)
tetrahedra = np.array([[0, 1, 2, 3]], dtype=np.int32)
mesh = mmgpy.Mesh(vertices, tetrahedra)
# From arrays - 2D triangular mesh
vertices_2d = np.array([
[0, 0],
[1, 0],
[0.5, 1],
], dtype=np.float64)
triangles = np.array([[0, 1, 2]], dtype=np.int32)
mesh_2d = mmgpy.Mesh(vertices_2d, triangles)
# From PyVista
import pyvista as pv
grid = pv.read("mesh.vtk")
mesh = mmgpy.Mesh(grid)
Checking Mesh Type¶
from mmgpy import MeshKind
mesh = mmgpy.Mesh(vertices, elements)
if mesh.kind == MeshKind.TETRAHEDRAL:
print("3D volume mesh")
elif mesh.kind == MeshKind.TRIANGULAR_2D:
print("2D planar mesh")
elif mesh.kind == MeshKind.TRIANGULAR_SURFACE:
print("3D surface mesh")
Accessing Mesh Data¶
# Get vertices and elements
vertices = mesh.get_vertices() # Shape: (n_vertices, 2 or 3)
triangles = mesh.get_triangles() # Shape: (n_triangles, 3)
tetrahedra = mesh.get_tetrahedra() # Shape: (n_tetrahedra, 4) - 3D only
# Get mesh statistics
size = mesh.get_mesh_size()
print(f"Vertices: {size['vertices']}")
Working with Fields¶
import numpy as np
# Set a scalar field (metric for sizing)
metric = np.ones(len(mesh.get_vertices()), dtype=np.float64) * 0.1
mesh["metric"] = metric
# Get a field
m = mesh["metric"]
Remeshing¶
from mmgpy import Mmg3DOptions
# With options object
opts = Mmg3DOptions(hmax=0.1, hausd=0.001)
result = mesh.remesh(opts)
# With keyword arguments
result = mesh.remesh(hmax=0.1, hausd=0.001)
# Convenience methods
result = mesh.remesh_optimize() # Quality only
result = mesh.remesh_uniform(size=0.1) # Uniform size
# Check results
print(f"Vertices: {result.vertices_before} -> {result.vertices_after}")
print(f"Quality: {result.quality_mean_before:.3f} -> {result.quality_mean_after:.3f}")
Local Sizing¶
# Set local refinement regions
mesh.set_size_sphere(center=[0.5, 0.5, 0.5], radius=0.2, size=0.01)
mesh.set_size_box(bounds=[[0, 0, 0], [0.3, 0.3, 0.3]], size=0.02)
# Remesh with local sizing applied
result = mesh.remesh(hmax=0.1)
# Clear sizing constraints
mesh.clear_local_sizing()
Validation¶
# Simple validation
if mesh.validate():
print("Mesh is valid")
# Detailed validation report
report = mesh.validate(detailed=True)
print(f"Valid: {report.is_valid}")
print(f"Quality: min={report.quality.min:.3f}, mean={report.quality.mean:.3f}")
# Strict mode (raises on issues)
mesh.validate(strict=True)