Skip to content

Mesh Classes

Deprecated in 0.12, removed in 0.13

mmgpy.Mesh and mmgpy.MeshCheckpoint are deprecated. New code should use the .mmg PyVista accessor (pv.read("foo.mesh").mmg.remesh(...)). See the migration guide for a method-by-method mapping. MeshKind survives in 0.13 and is still returned by dataset.mmg.kind.

Unified Mesh Class

mmgpy.Mesh

Mesh(
    source: NDArray[floating] | str | Path | UnstructuredGrid | PolyData,
    cells: NDArray[integer] | None = None,
    refs: NDArray[integer] | None = None,
    edges: NDArray[integer] | None = None,
    edge_refs: 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.

kind property

kind: MeshKind

Get the mesh kind.

Returns

MeshKind The type of mesh.

vtk property

vtk: UnstructuredGrid | PolyData

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() -> NDArray[np.float64]

Get vertex coordinates.

Returns

ndarray Vertex coordinates (Nx2 for 2D, Nx3 for 3D).

get_vertices_with_refs

get_vertices_with_refs() -> tuple[NDArray[np.float64], NDArray[np.int64]]

Get vertex coordinates and reference markers.

Returns

vertices : ndarray Vertex coordinates. refs : ndarray Reference markers for each vertex.

set_vertices

set_vertices(vertices: NDArray[float64], refs: NDArray[int64] | None = None) -> None

Set vertex coordinates.

Parameters

vertices : ndarray Vertex coordinates. refs : ndarray, optional Reference markers for each vertex.

get_triangles

get_triangles() -> NDArray[np.int32]

Get triangle connectivity.

Returns

ndarray Triangle connectivity (Nx3).

get_triangles_with_refs

get_triangles_with_refs() -> tuple[NDArray[np.int32], NDArray[np.int64]]

Get triangle connectivity and reference markers.

Returns

triangles : ndarray Triangle connectivity. refs : ndarray Reference markers for each triangle.

set_triangles

set_triangles(triangles: NDArray[int32], refs: NDArray[int64] | None = None) -> None

Set triangle connectivity.

Parameters

triangles : ndarray Triangle connectivity (Nx3). refs : ndarray, optional Reference markers for each triangle.

get_edges

get_edges() -> NDArray[np.int32]

Get edge connectivity.

Returns

ndarray Edge connectivity (Nx2).

get_edges_with_refs

get_edges_with_refs() -> tuple[NDArray[np.int32], NDArray[np.int64]]

Get edge connectivity and reference markers.

Returns

edges : ndarray Edge connectivity. refs : ndarray Reference markers for each edge.

set_edges

set_edges(edges: NDArray[int32], refs: NDArray[int64] | None = None) -> None

Set edge connectivity.

Parameters

edges : ndarray Edge connectivity (Nx2). refs : ndarray, optional Reference markers for each edge.

get_tetrahedra

get_tetrahedra() -> NDArray[np.int32]

Get tetrahedra connectivity.

Only available for TETRAHEDRAL meshes.

Returns

ndarray Tetrahedra connectivity (Nx4).

Raises

TypeError If mesh is not TETRAHEDRAL.

get_tetrahedra_with_refs

get_tetrahedra_with_refs() -> tuple[NDArray[np.int32], NDArray[np.int64]]

Get tetrahedra connectivity and reference markers.

Only available for TETRAHEDRAL meshes.

Returns

tetrahedra : ndarray Tetrahedra connectivity. refs : ndarray Reference markers for each tetrahedron.

Raises

TypeError If mesh is not TETRAHEDRAL.

get_elements

get_elements() -> NDArray[np.int32]

Get primary element connectivity (alias for get_tetrahedra).

Only available for TETRAHEDRAL meshes.

Returns

ndarray Element connectivity (Nx4 tetrahedra).

Raises

TypeError If mesh is not TETRAHEDRAL.

get_elements_with_refs

get_elements_with_refs() -> tuple[NDArray[np.int32], NDArray[np.int64]]

Get primary element connectivity and reference markers.

Only available for TETRAHEDRAL meshes.

Returns

elements : ndarray Element connectivity. refs : ndarray Reference markers for each element.

Raises

TypeError If mesh is not TETRAHEDRAL.

set_field

set_field(key: str, value: NDArray[float64]) -> None

Set a solution field.

.. deprecated:: 0.6.0 Use dictionary syntax instead: mesh[key] = value.

Parameters

key : str Field name. value : ndarray Field values (one per vertex).

get_field

get_field(key: str) -> NDArray[np.float64]

Get a solution field.

.. deprecated:: 0.6.0 Use dictionary syntax instead: value = mesh[key].

Parameters

key : str Field name.

Returns

ndarray Field values.

__setitem__

__setitem__(key: str, value: NDArray[float64]) -> None

Set a field using dictionary syntax.

Smart routing based on key and value shape:

  • mesh["metric"] = scalar → scalar metric (Nx1)
  • mesh["metric"] = tensor → anisotropic metric (Nx3 for 2D, Nx6 for 3D)
  • mesh["displacement"] = v → MMG displacement field
  • mesh["temperature"] = v → user-defined field (for transfer)
Parameters

key : str Field name. Known MMG fields ("metric", "displacement", "levelset", "tensor") are routed to the C++ bindings. All other names are stored as user fields. value : ndarray Field values (one per vertex).

Examples

mesh["metric"] = np.ones((n, 1)) * 0.1 # scalar metric mesh["metric"] = anisotropic_tensor # Nx6 tensor metric mesh["temperature"] = temp_array # user field

__getitem__

__getitem__(key: str) -> NDArray[np.float64]

Get a field using dictionary syntax.

Parameters

key : str Field name. Known MMG fields are retrieved from C++ bindings. All other names are retrieved from user fields.

Returns

ndarray Field values.

Raises

KeyError If the field does not exist or has not been set.

Examples

metric = mesh["metric"] temperature = mesh["temperature"]

__contains__

__contains__(key: str) -> bool

Check if a field exists.

Parameters

key : str Field name.

Returns

bool True if the field exists (either as an MMG field or user field).

Examples

"temperature" in mesh # True if user field was set "metric" in mesh # True if metric field is set

set_user_field

set_user_field(name: str, values: NDArray[float64]) -> None

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_field(name: str) -> NDArray[np.float64]

Get a user-defined field.

Parameters

name : str Field name.

Returns

ndarray Field values at vertices.

Raises

KeyError If the field does not exist.

get_user_fields

get_user_fields() -> dict[str, NDArray[np.float64]]

Get all user-defined fields.

Returns

dict[str, ndarray] Dictionary mapping field names to values.

clear_user_fields

clear_user_fields() -> None

Remove all user-defined fields.

has_user_field

has_user_field(name: str) -> bool

Check if a user field exists.

Parameters

name : str Field name.

Returns

bool True if the field exists.

get_bounds

get_bounds() -> tuple[NDArray[np.float64], NDArray[np.float64]]

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_center_of_mass() -> NDArray[np.float64]

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_volume() -> float

Compute the total volume of the mesh.

Only available for 3D volume meshes (TETRAHEDRAL).

Returns

float Total volume in mesh units cubed.

Raises

TypeError If mesh is not a 3D volume mesh.

Examples

mesh = Mesh(vertices, cells) volume = mesh.compute_volume() print(f"Volume: {volume:.2f} mm^3")

compute_surface_area

compute_surface_area() -> float

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_diagonal() -> float

Get the diagonal length of the bounding box.

Returns

float The diagonal length of the bounding box.

Raises

ValueError If the mesh has no vertices.

Examples

mesh = Mesh(vertices, cells) diagonal = mesh.get_diagonal() print(f"Bounding box diagonal: {diagonal:.2f}")

get_adjacent_elements

get_adjacent_elements(idx: int) -> NDArray[np.int32]

Get indices of elements adjacent to a given element.

Parameters

idx : int Element index (1-based for MMG).

Returns

ndarray Indices of adjacent elements.

get_vertex_neighbors

get_vertex_neighbors(idx: int) -> NDArray[np.int32]

Get indices of vertices connected to a given vertex.

Parameters

idx : int Vertex index (1-based for MMG).

Returns

ndarray Indices of neighboring vertices.

get_element_quality

get_element_quality(idx: int) -> float

Get quality metric for a single element.

Parameters

idx : int Element index (1-based for MMG).

Returns

float Quality metric (0-1, higher is better).

get_element_qualities

get_element_qualities() -> NDArray[np.float64]

Get quality metrics for all elements.

Returns

ndarray Quality metrics for all elements.

save

save(filename: str | Path) -> None

Save mesh to file.

Medit formats (.mesh, .meshb) are written directly by the MMG C library. All other formats (.vtk, .vtu, .vtp, .stl, .obj, ...) are converted via PyVista.

Parameters

filename : str or Path Output file path. Format determined by extension.

load_sol

load_sol(filename: str | Path) -> None

Load a Medit solution file (.sol/.solb) and set the field on the mesh.

Both text (.sol) and binary (.solb) formats are supported via the MMG C library.

Parameters

filename : str or Path Path to a .sol or .solb file.

save_sol

save_sol(filename: str | Path) -> None

Save the solution/metric field to a Medit .sol file.

Both text (.sol) and binary (.solb) formats are supported via the MMG C library.

Parameters

filename : str or Path Output path for the .sol or .solb file.

remesh

remesh(
    options: Mmg3DOptions | Mmg2DOptions | MmgSOptions | None = None,
    *,
    input_sol: str | Path | NDArray[float64] | 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. input_sol : str, Path, or ndarray, optional Solution data for adaptive remeshing. Can be a path to a Medit .sol file, or a numpy array (scalar metric, Nx3/Nx6 anisotropic tensor, or Nx2/Nx3 displacement vector). 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], **kwargs: Any) -> RemeshResult

Forward to :func:mmgpy.move_mesh (deprecated).

Historically this method called MMG's ELAS-bound Lagrangian path; the bundled MMG is built USE_ELAS=OFF, so it never produced anything useful. The shim now applies displacement and remeshes via the pure-Python Laplacian propagator (or the optional fedoo elasticity solver) and returns a minimal :class:RemeshResult. It will be removed in a future release; new code should call :func:mmgpy.move_mesh or dataset.mmg.move(...) directly.

Raises

TypeError If mesh is TRIANGULAR_SURFACE.

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

remesh_optimize(
    *, progress: ProgressParam = True, verbose: int | None = None
) -> RemeshResult

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_size_sphere(
    center: Sequence[float] | NDArray[float64], radius: float, size: float
) -> None

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_size_box(bounds: Sequence[Sequence[float]] | NDArray[float64], size: float) -> None

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.

clear_local_sizing

clear_local_sizing() -> None

Clear all local sizing constraints.

get_local_sizing_count

get_local_sizing_count() -> int

Get the number of local sizing constraints.

Returns

int Number of sizing constraints.

apply_local_sizing

apply_local_sizing() -> None

Apply local sizing constraints to the metric field.

This is called automatically before remeshing if sizing constraints have been added.

to_pyvista

to_pyvista(
    *, include_refs: bool = True, include_edges: bool = False
) -> pv.UnstructuredGrid | pv.PolyData

Convert to PyVista mesh.

Parameters

include_refs : bool Include reference markers as cell data. include_edges : bool Include MMG edges (ridges, boundary edges) as LINE cells in the output. Defaults to False so the returned mesh contains only the primary cell type (triangles or tetrahedra), which matches typical downstream code (matplotlib tripcolor, area computations, etc.). Set True for round-trip / file-save workflows that need to preserve edge markers.

Returns

pv.UnstructuredGrid | pv.PolyData PyVista mesh object.

plot

plot(*, show_edges: bool = True, **kwargs: Any) -> None

Plot the mesh using PyVista.

Parameters

show_edges : bool Show mesh edges (default: True). **kwargs : Any Additional arguments passed to PyVista's plot() method.

Examples

mesh = Mesh(vertices, cells) mesh.plot() # Simple plot with edges

mesh.plot(color="blue", opacity=0.8) # Custom styling

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__

__enter__() -> Mesh

Enter the context manager.

Returns

Mesh The mesh instance.

Examples

with Mesh(vertices, cells) as mesh: ... mesh.remesh(hmax=0.1) ... mesh.save("output.vtk")

__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

checkpoint() -> MeshCheckpoint

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

copy() -> Generator[Mesh, None, None]

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_from(other: Mesh) -> None

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

Legacy Usage Examples

The snippets below illustrate the deprecated Mesh API. For new code, see pyvista-integration.md and the migration guide.

Creating Meshes

import mmgpy
import numpy as np

mesh = mmgpy.Mesh("input.mesh")  # deprecated — use pv.read("input.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)  # deprecated

Checking Mesh Type

from mmgpy import MeshKind

mesh = mmgpy.Mesh(vertices, tetrahedra)  # deprecated

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")

Remeshing

from mmgpy import Mmg3DOptions

opts = Mmg3DOptions(hmax=0.1, hausd=0.001)
result = mesh.remesh(opts)  # deprecated; prefer pv_dataset.mmg.remesh(opts)

Validation

report = mesh.validate(detailed=True)  # deprecated; prefer dataset.mmg.validate(detailed=True)