Skip to content

I/O Functions

This page documents the input/output functions for loading and saving meshes.

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 meshio (40+ formats including VTK, Gmsh, STL, OBJ, etc.) 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 that meshio does not understand.

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 (native MMG loading, preserves all MMG keywords)
  • VTK: .vtk, .vtu, .vtp (via meshio)
  • Gmsh: .msh (via meshio)
  • STL: .stl (via meshio)
  • OBJ: .obj (via meshio)
  • PLY: .ply (via meshio)
  • And many more via meshio...
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

Supported Formats

mmgpy supports 40+ file formats via meshio:

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 mmgpy

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

# Returns unified Mesh object
print(f"Type: {mesh.kind}")  # MeshKind.TETRAHEDRAL, etc.

PyVista Conversion

mmgpy.from_pyvista

from_pyvista(
    mesh: UnstructuredGrid | PolyData, mesh_type: type[MmgMesh3D]
) -> MmgMesh3D
from_pyvista(
    mesh: UnstructuredGrid | PolyData, mesh_type: type[MmgMesh2D]
) -> MmgMesh2D
from_pyvista(mesh: UnstructuredGrid | PolyData, mesh_type: type[MmgMeshS]) -> MmgMeshS
from_pyvista(
    mesh: UnstructuredGrid | PolyData, mesh_type: None = None
) -> MmgMesh3D | MmgMesh2D | MmgMeshS
from_pyvista(
    mesh: UnstructuredGrid | PolyData,
    mesh_type: type[MmgMesh3D | MmgMesh2D | MmgMeshS] | None = None,
) -> MmgMesh3D | MmgMesh2D | MmgMeshS

Convert a PyVista mesh to an mmgpy mesh.

Parameters:

  • mesh (UnstructuredGrid | PolyData) –

    PyVista mesh (UnstructuredGrid or PolyData).

  • mesh_type (type[MmgMesh3D | MmgMesh2D | MmgMeshS] | None, default: None ) –

    Target mesh class. If None, auto-detects based on: - UnstructuredGrid with tetrahedra → MmgMesh3D - PolyData with 2D points (z~=0) -> MmgMesh2D - PolyData with 3D points → MmgMeshS

Returns:

  • MmgMesh3D | MmgMesh2D | MmgMeshS

    The appropriate mmgpy mesh instance.

Raises:

  • ValueError

    If mesh type cannot be determined or is incompatible.

Note

When auto-detecting mesh type for PolyData, a mesh is considered 2D (and converted to MmgMesh2D) if all z-coordinates are within 1e-8 of zero. For thin 3D meshes near z=0, explicitly specify mesh_type=MmgMeshS.

Example

import pyvista as pv from mmgpy import from_pyvista, MmgMeshS

Auto-detect mesh type

grid = pv.read("tetra_mesh.vtk") mesh3d = from_pyvista(grid)

Explicit mesh type for thin 3D surfaces

surface = pv.read("surface.stl") mesh_s = from_pyvista(surface, MmgMeshS)

options: show_root_heading: true

mmgpy.to_pyvista

to_pyvista(mesh: MmgMesh3D, *, include_refs: bool = True) -> pv.UnstructuredGrid
to_pyvista(mesh: MmgMesh2D, *, include_refs: bool = True) -> pv.PolyData
to_pyvista(mesh: MmgMeshS, *, include_refs: bool = True) -> pv.PolyData
to_pyvista(
    mesh: MmgMesh3D | MmgMesh2D | MmgMeshS, *, include_refs: bool = True
) -> pv.UnstructuredGrid | pv.PolyData

Convert an mmgpy mesh to a PyVista mesh.

Parameters:

  • mesh (MmgMesh3D | MmgMesh2D | MmgMeshS) –

    mmgpy mesh instance (MmgMesh3D, MmgMesh2D, or MmgMeshS).

  • include_refs (bool, default: True ) –

    If True, include element references as cell_data.

Returns:

  • UnstructuredGrid | PolyData

    PyVista mesh: - MmgMesh3D → UnstructuredGrid with tetrahedra - MmgMesh2D → PolyData with triangular faces (z=0) - MmgMeshS → PolyData with triangular faces

Raises:

  • TypeError

    If mesh is not an mmgpy mesh type.

Example

from mmgpy import MmgMesh3D, to_pyvista

mesh = MmgMesh3D(vertices, elements) mesh.remesh(hmax=0.1) grid = to_pyvista(mesh) grid.plot()

options: show_root_heading: true

From PyVista

import mmgpy
import pyvista as pv

# Create PyVista geometry
sphere = pv.Sphere(radius=1.0)

# Convert to surface mesh
mesh = mmgpy.from_pyvista(sphere, mesh_type="surface")

# For volumetric meshes (needs tetrahedral cells)
volume = pv.Box().triangulate().delaunay_3d()
mesh_3d = mmgpy.from_pyvista(volume, mesh_type="3d")

# For 2D meshes
plane = pv.Plane()
mesh_2d = mmgpy.from_pyvista(plane, mesh_type="2d")

To PyVista

import mmgpy

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

# Convert to PyVista
pv_mesh = mesh.to_pyvista()

# Or using function
pv_mesh = mmgpy.to_pyvista(mesh)

# Visualize
pv_mesh.plot(show_edges=True)

Saving Meshes

Meshes are saved using the save() method:

mesh.save("output.mesh")  # MMG native
mesh.save("output.vtk")   # VTK format
mesh.save("output.stl")   # STL (surface only)

Format is inferred from the file extension.

Direct Class Loading

Each mesh class can also load from files directly:

import mmgpy

# 3D mesh
mesh_3d = mmgpy.Mesh("volume.mesh")

# 2D mesh
mesh_2d = mmgpy.Mesh("planar.mesh")

# Surface mesh
mesh_s = mmgpy.Mesh("surface.stl")

Complete Example

import mmgpy
import pyvista as pv

# Load from file
mesh = mmgpy.read("input.mesh")
print(f"Loaded {mesh.kind} mesh")

# Remesh
mesh.remesh(hmax=0.1)

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

# Convert to PyVista for visualization
pv_mesh = mesh.to_pyvista()
pv_mesh.save("output_pv.vtk")

# Or create from PyVista
torus = pv.ParametricTorus()
torus_mesh = mmgpy.from_pyvista(torus, mesh_type="surface")
torus_mesh.remesh(hmax=0.1)
torus_mesh.save("torus.mesh")

Tips

  1. MMG native format: Use .mesh for fastest I/O with MMG
  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