lapy.heat

Functions for computing heat kernel and diffusion.

Inputs are eigenvalues and eigenvectors (for heat kernel) and the mesh geometries (tet or tria mesh) for heat diffusion.

lapy.heat.diagonal(t, x, evecs, evals, n)[source]

Compute heat kernel diagonal ( K(t,x,x,) ).

For a given time t (can be a vector) using only the first n smallest eigenvalues and eigenvectors.

Parameters:
tfloat | array

Time or a row vector of time values.

xarray

Vertex ids for the positions of K(t,x,x).

evecsarray

Eigenvectors (matrix: vnum x evecsnum).

evalsarray

Vector of eigenvalues (col vector: evecsnum x 1).

nint

Number of evecs and vals to use (smaller or equal length).

Returns:
harray

Matrix, rows: vertices selected in x, cols: times in t.

lapy.heat.diffusion(geometry, vids, m=1.0, aniso=None, use_cholmod=False)[source]

Compute the heat diffusion from initial vertices in vids.

It uses the backward Euler solution \(t = m l^2\), where l describes the average edge length.

Parameters:
geometryTriaMesh | TetMesh

Geometric object on which to run diffusion.

vidsarray

Vertex index or indices where initial heat is applied.

mfloat, default=1.0

Factor to compute time of heat evolution.

anisoint

Number of smoothing iterations for curvature computation on vertices.

use_cholmodbool, default=False
Which solver to use:
  • True : Use Cholesky decomposition from scikit-sparse cholmod.

  • False: Use spsolve (LU decomposition).

Returns:
vfunc: array of shape (n, 1)

Heat diffusion at vertices.

lapy.heat.kernel(t, vfix, evecs, evals, n)[source]

Compute heat kernel from all points to a fixed point (vfix).

For a given time t (using only the first n smallest eigenvalues and eigenvectors):

\[K_t (p,q) = \sum_j \ exp(-eval_j \ t) \ evec_j(p) \ evec_j(q)\]
Parameters:
tfloat | array

Time (can also be a row vector, if passing multiple times).

vfixarray

Fixed vertex index.

evecsarray

Matrix of eigenvectors (M x N), M=#vertices, N=#eigenvectors.

evalsarray

Column vector of eigenvalues (N).

nint

Number of eigenvalues/vectors used in heat kernel (n<=N).

Returns:
harray

Matrix m rows: all vertices, cols: times in t.