Skip to content

I/O

This page documents how to load and save meshes with mmgpy.

mmgpy registers a Medit reader/writer plugin with PyVista on import, so pv.read("foo.mesh") and dataset.save("foo.mesh") both go through MMG's native I/O. Other formats are handled by PyVista directly (which uses meshio under the hood).

Reading Meshes

mmgpy.read

read(
    source: str | Path | UnstructuredGrid | PolyData, mesh_kind: MeshKind | None = None
) -> Mesh

Read a mesh from a file or PyVista object.

This function provides unified mesh loading from any format supported by PyVista or directly from PyVista mesh objects.

For Medit format (.mesh) files, native MMG loading is used to preserve MMG-specific keywords like Ridges, RequiredVertices, Tangents, and reference markers.

Parameters:

  • source (str | Path | UnstructuredGrid | PolyData) –

    File path (str or Path) or PyVista mesh object.

  • mesh_kind (MeshKind | None, default: None ) –

    Force a specific mesh kind instead of auto-detection. - MeshKind.TETRAHEDRAL: 3D volumetric mesh - MeshKind.TRIANGULAR_2D: 2D planar mesh - MeshKind.TRIANGULAR_SURFACE: 3D surface mesh - None: Auto-detect based on element types and coordinates

Returns:

  • Mesh

    A Mesh instance with the appropriate kind.

Raises:

  • ValueError

    If mesh kind cannot be determined or file cannot be read.

  • TypeError

    If source type is not supported.

  • FileNotFoundError

    If file does not exist.

Auto-detection logic
  • Has tetrahedra → TETRAHEDRAL
  • Has triangles + 3D coords → TRIANGULAR_SURFACE
  • Has triangles + 2D coords (or z~=0) -> TRIANGULAR_2D
Supported file formats
  • Medit: .mesh, .meshb (native MMG loading, preserves all MMG keywords)
  • VTK: .vtk, .vtu, .vtp (via PyVista)
  • STL: .stl (via PyVista)
  • OBJ: .obj (via PyVista)
  • PLY: .ply (via PyVista)
  • Gmsh, Abaqus, NASTRAN, Exodus, XDMF, etc. (via PyVista + meshio)
  • And many more via PyVista (which uses meshio as fallback)...
Example

import mmgpy

Auto-detect mesh kind from file

mesh = mmgpy.read("tetra_mesh.vtk") mesh.kind # MeshKind.TETRAHEDRAL

Force specific mesh kind

mesh = mmgpy.read("mesh.vtk", mesh_kind=MeshKind.TRIANGULAR_SURFACE)

Read Medit file with native loading (preserves Ridges, etc.)

mesh = mmgpy.read("mesh.mesh")

Read from PyVista object

import pyvista as pv grid = pv.read("mesh.vtk") mesh = mmgpy.read(grid)

options: show_root_heading: true

mmgpy.read returns a Mesh (deprecated). For new code, use pv.read(...) and the .mmg accessor.

Supported Formats

Format Extensions Notes
MMG native .mesh, .meshb Recommended for MMG
VTK Legacy .vtk Universal, ParaView compatible
VTK XML .vtu, .vtp Modern VTK format
STL .stl Surface meshes only
OBJ .obj Surface meshes only
PLY .ply Point cloud / mesh
GMSH .msh Popular for FEM
Abaqus .inp FEM format
CGNS .cgns CFD format
Exodus II .e, .exo Sandia format
ANSYS .ansys FEM format
MED .med Salome format
And many more... See meshio documentation

Usage

import pyvista as pv
import mmgpy  # noqa: F401  -- registers reader/writer + accessor

# Auto-detect format from extension
mesh = pv.read("input.mesh")
mesh = pv.read("input.vtk")
mesh = pv.read("input.stl")

print(f"Kind: {mesh.mmg.kind}")  # MeshKind.TETRAHEDRAL, etc.

From PyVista Primitives

The accessor works on any pv.UnstructuredGrid or pv.PolyData, so you don't need an explicit conversion:

import pyvista as pv
import mmgpy  # noqa: F401

# Surface
sphere = pv.Sphere(radius=1.0)
remeshed_surface = sphere.mmg.remesh(hsiz=0.1)

# Volume (needs tetrahedral cells)
volume = pv.Box().triangulate().delaunay_3d()
remeshed_volume = volume.mmg.remesh(hsiz=0.2)

# 2D plane
plane = pv.Plane()
remeshed_plane = plane.mmg.remesh(hsiz=0.1)

Saving Meshes

Use PyVista's save():

# MMG native (registered reader/writer)
remeshed_volume.save("output.mesh")

# VTK formats (handled by PyVista)
remeshed_volume.save("output.vtk")
remeshed_volume.save("output.vtu")

# Surface meshes also support .stl
remeshed_surface.save("output.stl")

Format is inferred from the file extension.

Complete Example

import pyvista as pv
import mmgpy  # noqa: F401

mesh = pv.read("input.mesh")
print(f"Loaded {mesh.mmg.kind} mesh")

remeshed = mesh.mmg.remesh(hmax=0.1)

# Save to different formats
remeshed.save("output.mesh")   # MMG native (fast)
remeshed.save("output.vtk")    # For ParaView
remeshed.save("output.vtu")    # VTK XML format

# Or build from a PyVista primitive
torus = pv.ParametricTorus()
remeshed_torus = torus.mmg.remesh(hmax=0.1)
remeshed_torus.save("torus.mesh")

Tips

  1. MMG native format: Use .mesh for fastest I/O with MMG (mmgpy's reader plugin handles ridges and reference markers natively).
  2. VTK for visualization: Use .vtk or .vtu for ParaView.
  3. Surface formats: STL and OBJ are surface-only.
  4. Binary formats: Some formats support binary (faster, smaller).
  5. Field data: Most formats preserve scalar/vector fields via point_data/cell_data.