lapy.shapedna¶
Functions for computing and comparing Laplace spectra.
Includes code for solving the anisotropic Laplace-Beltrami eigenvalue problem as well as functions for normalization and comparison of Laplace spectra (shapeDNA descriptors).
The shapeDNA is a descriptor based on the eigenvalues and eigenvectors of the Laplace-Beltrami operator and can be used for shape analysis and comparison.
- lapy.shapedna.compute_distance(ev1, ev2, dist='euc')[source]¶
Compute the shape dissimilarity from two shapeDNA descriptors.
Computes a distance metric between two sets of eigenvalues to quantify the dissimilarity between two shapes.
- Parameters:
- ev1
np.ndarray First set of sorted eigenvalues, shape (k,).
- ev2
np.ndarray Second set of sorted eigenvalues, shape (k,).
- dist{‘euc’}, default=’euc’
Distance measure. Currently only ‘euc’ (Euclidean) is implemented.
- ev1
- Returns:
floatDistance between the eigenvalue arrays.
- Raises:
ValueErrorIf dist is not ‘euc’ (other distance metrics not yet implemented).
- Parameters:
- Return type:
Notes
The eigenvalue arrays should have the same length and be normalized and reweighted in the same way for meaningful comparison.
- lapy.shapedna.compute_shapedna(geom, k=50, lump=False, aniso=None, aniso_smooth=10, use_cholmod=False)[source]¶
Compute the shapeDNA descriptor for triangle or tetrahedral meshes.
The shapeDNA descriptor consists of the eigenvalues and eigenvectors of the Laplace-Beltrami operator and can be used for shape analysis and comparison.
- Parameters:
- geom
TriaMesh|TetMesh Mesh geometry.
- k
int, default=50 Number of eigenvalues/eigenvectors to compute.
- lump
bool, default=False If True, lump the mass matrix (diagonal). See
lapy.Solverclass.- aniso
float|tupleof shape (2,) |None, default=None Anisotropy for curvature-based anisotropic Laplace. See
lapy.Solverclass.- aniso_smooth
int, default=10 Number of smoothing iterations for curvature computation on vertices. See
lapy.Solverclass.- use_cholmod
bool, default=False If True, attempts to use the Cholesky decomposition for improved execution speed. Requires the
scikit-sparselibrary. If it cannot be found, an error will be thrown. If False, will use slower LU decomposition.
- geom
- Returns:
dictA dictionary with the following keys:
‘Refine’ : int - Refinement level (0)
‘Degree’ : int - Polynomial degree (1)
‘Dimension’ : int - Mesh dimension (2 for TriaMesh, 3 for TetMesh)
‘Elements’ : int - Number of mesh elements
‘DoF’ : int - Degrees of freedom (number of vertices)
‘NumEW’ : int - Number of eigenvalues computed
‘Eigenvalues’ : np.ndarray - Array of eigenvalues, shape (k,)
‘Eigenvectors’ : np.ndarray - Array of eigenvectors, shape (n_vertices, k)
- Parameters:
- Return type:
- lapy.shapedna.normalize_ev(geom, evals, method='geometry')[source]¶
Normalize eigenvalues for a 2D surface or a 3D solid.
Normalizes eigenvalues to unit surface area or unit volume, enabling meaningful comparison between different shapes.
- Parameters:
- geom
TriaMeshorTetMesh Mesh geometry.
- evals
np.ndarray Set of sorted eigenvalues, shape (k,).
- method{‘surface’, ‘volume’, ‘geometry’}, default=’geometry’
Normalization method:
‘surface’: Normalize to unit surface area
‘volume’: Normalize to unit volume
‘geometry’: Automatically choose surface for TriaMesh, volume for TetMesh
- geom
- Returns:
np.ndarrayVector of re-weighted eigenvalues, shape (k,).
- Raises:
ValueErrorIf the method is not one of ‘surface’, ‘volume’, or ‘geometry’. If the geometry type is unsupported for the chosen normalization. If method=volume and the volume of a surface is not defined.
- Parameters:
- Return type:
Notes
For TriaMesh with ‘volume’ method, the mesh must be closed and oriented to compute a valid enclosed volume.
- lapy.shapedna.reweight_ev(evals)[source]¶
Apply linear re-weighting to eigenvalues.
Divides each eigenvalue by its index to reduce the influence of higher eigenvalues, which tend to be less stable.
- Parameters:
- evals
np.ndarray Set of sorted eigenvalues, shape (k,).
- evals
- Returns:
np.ndarrayVector of re-weighted eigenvalues, shape (k,). Each eigenvalue is divided by its 1-based index:
evals[i] / (i+1).
- Parameters:
evals (ndarray)
- Return type:
Notes
This reweighting scheme gives less importance to higher eigenvalues, which are typically more sensitive to discretization and numerical errors.