Skip to content

Metrics

This page documents the metric tensor operations in the mmgpy.metrics module.

Overview

Metric tensors control anisotropic mesh adaptation. A metric at each vertex specifies:

  • Isotropic: Target edge length (single scalar)
  • Anisotropic: Target lengths along different directions (tensor)

Metric Creation

mmgpy.metrics.create_isotropic_metric

create_isotropic_metric(
    h: float | NDArray[float64], n_vertices: int | None = None, dim: int = 3
) -> NDArray[np.float64]

Create an isotropic metric field from scalar sizing values.

Parameters

h : float or array_like Desired element size(s). If scalar, same size at all vertices. If array, must have shape (n_vertices,) or (n_vertices, 1). n_vertices : int, optional Number of vertices. Required if h is a scalar. dim : int, optional Mesh dimension (2 or 3). Default is 3.

Returns

NDArray[np.float64] Metric tensor array with shape (n_vertices, n_components) where n_components is 6 for 3D and 3 for 2D.

Examples

metric = create_isotropic_metric(0.1, n_vertices=100, dim=3) metric.shape (100, 6)

sizes = np.linspace(0.1, 0.5, 100) metric = create_isotropic_metric(sizes, dim=3) metric.shape (100, 6)

options: show_root_heading: true

mmgpy.metrics.create_anisotropic_metric

create_anisotropic_metric(
    sizes: NDArray[float64], directions: NDArray[float64] | None = None
) -> NDArray[np.float64]

Create an anisotropic metric tensor from principal sizes and directions.

The metric tensor M is constructed as: M = R @ D @ R.T where D is diagonal with D[i,i] = 1/sizes[i]^2 and R contains the principal direction vectors as columns.

Parameters

sizes : array_like Principal element sizes. Shape (3,) for 3D, (2,) for 2D. Can also be (n_vertices, 3) or (n_vertices, 2) for per-vertex sizes. directions : array_like, optional Principal direction vectors. Shape (3, 3) or (2, 2) for single metric, or (n_vertices, 3, 3) or (n_vertices, 2, 2) for per-vertex directions. Columns are eigenvectors. If None, uses identity (coordinate-aligned).

Returns

NDArray[np.float64] Metric tensor(s). Shape (6,) for single 3D metric, (3,) for single 2D, or (n_vertices, 6) / (n_vertices, 3) for per-vertex metrics.

Examples

Create a metric with 10x stretch in x-direction:

sizes = np.array([0.1, 1.0, 1.0]) # Small in x, large in y,z metric = create_anisotropic_metric(sizes) metric array([100., 0., 0., 1., 0., 1.])

Create a rotated metric:

import numpy as np theta = np.pi / 4 # 45 degrees R = np.array([[np.cos(theta), -np.sin(theta), 0], ... [np.sin(theta), np.cos(theta), 0], ... [0, 0, 1]]) sizes = np.array([0.1, 1.0, 1.0]) metric = create_anisotropic_metric(sizes, R)

options: show_root_heading: true

mmgpy.metrics.create_metric_from_hessian

create_metric_from_hessian(
    hessian: NDArray[float64],
    target_error: float = 0.001,
    hmin: float | None = None,
    hmax: float | None = None,
) -> NDArray[np.float64]

Create metric tensor from Hessian matrix for interpolation error control.

Given a Hessian H of a solution field, constructs a metric M such that the interpolation error is bounded by target_error. This is used for solution-adaptive mesh refinement.

The metric eigenvalues are: lambda_i = c * |hessian_eigenvalue_i| / target_error where c is a constant depending on the interpolation order.

Parameters

hessian : array_like Hessian tensor(s). Shape (6,) or (n, 6) for 3D, (3,) or (n, 3) for 2D. Components: [H11, H12, H13, H22, H23, H33] for 3D. target_error : float, optional Target interpolation error. Default is 1e-3. hmin : float, optional Minimum element size. Limits maximum metric eigenvalues. hmax : float, optional Maximum element size. Limits minimum metric eigenvalues.

Returns

NDArray[np.float64] Metric tensor(s) for adaptive remeshing.

Notes

For P1 interpolation, the interpolation error is bounded by: e <= (1/8) * h^2 * |d²u/ds²|_max

This function computes the metric that achieves a specified error bound.

options: show_root_heading: true

Metric Operations

mmgpy.metrics.intersect_metrics

intersect_metrics(
    m1: NDArray[float64], m2: NDArray[float64], dim: int | None = None
) -> NDArray[np.float64]

Compute the intersection of two metric tensors.

The intersection produces a metric that is at least as refined as both input metrics in all directions. This is useful for combining metrics from different sources (e.g., boundary layer + feature-based).

The intersection is computed via simultaneous diagonalization: M_intersect = M1^(1/2) @ N @ M1^(1/2) where N is diagonal with max eigenvalues of M1^(-1/2) @ M2 @ M1^(-1/2).

Parameters

m1, m2 : array_like Metric tensors to intersect. Must have same shape. Shape (6,) or (n, 6) for 3D, (3,) or (n, 3) for 2D. dim : int, optional Dimension (2 or 3). Inferred from tensor shape if not provided.

Returns

NDArray[np.float64] Intersected metric tensor(s), same shape as inputs.

options: show_root_heading: true

mmgpy.metrics.compute_metric_eigenpairs

compute_metric_eigenpairs(
    tensor: NDArray[float64], dim: int | None = None
) -> tuple[NDArray[np.float64], NDArray[np.float64]]

Extract principal sizes and directions from metric tensor(s).

Parameters

tensor : array_like Metric tensor(s). Shape (6,) or (n, 6) for 3D, (3,) or (n, 3) for 2D. dim : int, optional Dimension (2 or 3). Inferred from tensor shape if not provided.

Returns

tuple[NDArray, NDArray] (sizes, directions) where: - sizes: Principal element sizes, shape (3,) or (n, 3) for 3D - directions: Eigenvector matrices, shape (3, 3) or (n, 3, 3) for 3D Columns are eigenvectors corresponding to sizes.

Examples

tensor = np.array([100., 0., 0., 1., 0., 1.]) # 10x stretch in x sizes, directions = compute_metric_eigenpairs(tensor) sizes array([0.1, 1. , 1. ])

options: show_root_heading: true

Tensor Utilities

mmgpy.metrics.tensor_to_matrix

tensor_to_matrix(
    tensor: NDArray[float64], dim: int | None = None
) -> NDArray[np.float64]

Convert tensor storage format to full symmetric matrix.

Parameters

tensor : array_like Tensor in storage format. Shape (6,) or (n, 6) for 3D, (3,) or (n, 3) for 2D. dim : int, optional Dimension (2 or 3). Inferred from tensor shape if not provided.

Returns

NDArray[np.float64] Full symmetric matrix. Shape (3, 3) or (n, 3, 3) for 3D, (2, 2) or (n, 2, 2) for 2D.

options: show_root_heading: true

mmgpy.metrics.matrix_to_tensor

matrix_to_tensor(M: NDArray[float64]) -> NDArray[np.float64]

Convert full symmetric matrix to tensor storage format.

Parameters

M : array_like Full symmetric matrix. Shape (3, 3) or (n, 3, 3) for 3D, (2, 2) or (n, 2, 2) for 2D.

Returns

NDArray[np.float64] Tensor in storage format. Shape (6,) or (n, 6) for 3D, (3,) or (n, 3) for 2D.

options: show_root_heading: true

mmgpy.metrics.validate_metric_tensor

validate_metric_tensor(
    tensor: NDArray[float64], dim: int | None = None, *, raise_on_invalid: bool = True
) -> tuple[bool, str]

Validate that metric tensor(s) are positive-definite.

A valid metric tensor must be symmetric positive-definite, meaning all eigenvalues must be strictly positive.

Parameters

tensor : array_like Tensor(s) to validate. Shape (6,) or (n, 6) for 3D, (3,) or (n, 3) for 2D. dim : int, optional Dimension (2 or 3). Inferred from tensor shape if not provided. raise_on_invalid : bool, optional If True, raises ValueError on invalid tensors. Default is True.

Returns

tuple[bool, str] (is_valid, message) tuple.

Raises

ValueError If raise_on_invalid is True and tensor is not valid.

Examples

valid_tensor = np.array([1.0, 0.0, 0.0, 1.0, 0.0, 1.0]) validate_metric_tensor(valid_tensor) (True, 'Valid positive-definite metric tensor')

invalid_tensor = np.array([-1.0, 0.0, 0.0, 1.0, 0.0, 1.0]) validate_metric_tensor(invalid_tensor, raise_on_invalid=False) (False, 'Tensor has non-positive eigenvalues...')

options: show_root_heading: true

Usage Examples

Isotropic Metric

Create a metric for uniform element sizes:

import mmgpy
import mmgpy.metrics as metrics
import numpy as np

mesh = mmgpy.read("input.mesh")
n_vertices = mesh.get_mesh_size()["vertices"]

# Uniform size everywhere
sizes = np.ones(n_vertices) * 0.1
metric = metrics.create_isotropic_metric(sizes)

# Apply to mesh
mesh["metric"] = metric

# Remesh using the metric
result = mesh.remesh()

Variable Size Metric

Size varying with position:

import numpy as np

vertices = mesh.get_vertices()

# Size increases with distance from origin
distances = np.linalg.norm(vertices, axis=1)
sizes = 0.01 + 0.1 * distances

metric = metrics.create_isotropic_metric(sizes)
mesh["metric"] = metric

Anisotropic Metric

Different sizes in different directions:

import numpy as np

n_vertices = mesh.get_mesh_size()["vertices"]

# Define principal directions and sizes at each vertex
# directions: (n_vertices, 3, 3) - orthonormal basis at each vertex
# sizes: (n_vertices, 3) - sizes along each principal direction

# Example: stretch along z-axis
directions = np.tile(np.eye(3), (n_vertices, 1, 1))  # Identity basis
sizes = np.tile([0.1, 0.1, 0.05], (n_vertices, 1))   # Smaller in z

metric = metrics.create_anisotropic_metric(directions, sizes)
mesh["metric"] = metric

Metric from Hessian

Adapt mesh to solution curvature:

# Solution field (e.g., temperature)
solution = np.sin(vertices[:, 0] * 2 * np.pi)

# Compute Hessian (second derivatives) - requires additional computation
# hessian shape: (n_vertices, 6) for symmetric 3x3 tensor
hessian = compute_hessian(solution, mesh)  # Implementation needed

# Create metric from Hessian
metric = metrics.create_metric_from_hessian(
    hessian,
    epsilon=0.01,  # Target interpolation error
)

mesh["metric"] = metric

Metric Intersection

Combine multiple metrics (minimum size wins):

# Two different metrics
metric1 = metrics.create_isotropic_metric(sizes1)
metric2 = metrics.create_isotropic_metric(sizes2)

# Intersect: take minimum size in all directions
combined = metrics.intersect_metrics(metric1, metric2)
mesh["metric"] = combined

Extracting Metric Information

# Get current metric
metric = mesh["metric"]

# Extract eigenvalues and eigenvectors
eigenvalues, eigenvectors = metrics.compute_metric_eigenpairs(metric)

# eigenvalues shape: (n_vertices, 3) - 1/size^2 along each direction
# eigenvectors shape: (n_vertices, 3, 3) - principal directions

# Convert to sizes
sizes = 1.0 / np.sqrt(eigenvalues)
print(f"Size range: {sizes.min():.4f} to {sizes.max():.4f}")

Tensor Format Conversion

MMG uses symmetric tensors in Voigt notation:

# 3D: 6 components per vertex
# [M11, M12, M13, M22, M23, M33]

# 2D: 3 components per vertex
# [M11, M12, M22]

# Convert tensor to full matrix
tensor = mesh["metric"][0]  # First vertex
matrix = metrics.tensor_to_matrix(tensor)
print(matrix.shape)  # (3, 3)

# Convert matrix back to tensor
tensor_back = metrics.matrix_to_tensor(matrix)

Validation

Check metric tensor validity:

metric = mesh["metric"]

# Validate: must be symmetric positive definite
is_valid = metrics.validate_metric_tensor(metric)
if not is_valid:
    print("Warning: invalid metric tensor")

Metric Formats

3D Metrics (TETRAHEDRAL)

Symmetric 3x3 tensor stored as 6 components:

    [M11  M12  M13]
M = [M12  M22  M23]  -> [M11, M12, M13, M22, M23, M33]
    [M13  M23  M33]

Metric field shape: (n_vertices, 6)

2D Metrics (TRIANGULAR_2D)

Symmetric 2x2 tensor stored as 3 components:

    [M11  M12]
M = [M12  M22]  -> [M11, M12, M22]

Metric field shape: (n_vertices, 3)

Surface Metrics (TRIANGULAR_SURFACE)

Same as 3D: (n_vertices, 6)

Tips

  1. Isotropic first: Start with isotropic metrics, add anisotropy only when needed

  2. Size bounds: Ensure metric sizes are within reasonable bounds relative to domain size

  3. Gradation: MMG's hgrad parameter controls size gradation regardless of metric

  4. Validation: Always validate metric tensors before remeshing

  5. Combination: Use intersect_metrics to combine sizing from different sources