Matrix operators
Connectivity operators
adjacency_matrix(mesh, weights='one')
Computes and returns the adjacency matrix of a mesh, that is a (sparse) matrix M such that
M[i,j] = M[j,i] = weights[i,j] if (i,j) is an edge of the mesh
M[i,j] = 0 otherwise
Parameters:
Name | Type | Description | Default |
---|---|---|---|
mesh
|
Mesh
|
the input mesh |
required |
weights
|
str
|
How to weight the edges of the matrix. Options are: - "one" : every edge is 1 - "length" : every edge has weight corresponding to its length - a custom dict edge_id:int -> weight:float Defaults to "one". |
'one'
|
Returns:
Type | Description |
---|---|
coo_matrix
|
scipy.sparse.coo_matrix |
vertex_to_edge_operator(mesh, oriented=False)
Vertices to edges operator. Matrix M of size |V|x|E| where:
M[v,e] = 1 if and only if v is one extremity of edge e.
if oriented is set, M[v,e] = -1 if v is the origin of the edge and 1 if it is the arrival.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
mesh
|
Mesh
|
the input mesh |
required |
oriented
|
bool
|
whether to consider oriented edge and signed values or not. Defaults to False. |
False
|
Returns:
Type | Description |
---|---|
csc_matrix
|
scipy.sparse.csc_matrix |
vertex_to_face_operator(mesh)
Vertices to face operator. Matrix M of size |V| x |F| where:
M[v,f] = 1/len(f) if and only if v is one of the vertices of f
Parameters:
Name | Type | Description | Default |
---|---|---|---|
mesh
|
SurfaceMesh
|
input surface mesh |
required |
Returns:
Type | Description |
---|---|
csc_matrix
|
scipy.sparse.csc_matrix |
See also
mouette.attributes.interpolate_vertices_to_faces
Gradient operator
G = M.operators.gradient(mesh)
my_fun = mesh.vertices.get_attribute("f").as_array()
grad = G @ my_fun
gradient(mesh, conn=None, as_complex=True)
Computes the gradient operator, i.e. a |F| x |V| matrix G such that for any scalar function f defined over vertices of a surface mesh, Gf is its gradient inside faces.
Gf maps to each face of the mesh either a vector of \(\mathbb{R}^2\) or a complex number representing the gradient vector inside this face in local base.
See https://github.com/GCoiffier/mouette/blob/main/examples/gradient.py
Parameters:
Name | Type | Description | Default |
---|---|---|---|
mesh
|
SurfaceMesh
|
The input mesh |
required |
conn
|
SurfaceConnectionFaces
|
Connection objects specifying local bases of faces. Will be computed if not provided. Defaults to None. |
None
|
as_complex
|
bool
|
whether the output is |F| complex values of 2|F| float values |
True
|
Raises:
Type | Description |
---|---|
Exception
|
Fails if the mesh is not a triangulation |
Returns:
Type | Description |
---|---|
csc_matrix
|
scipy.sparse.csc_matrix: Gradient operator |
Laplacian operators
https://en.wikipedia.org/wiki/Discrete_Laplace_operator#Mesh_Laplacians
We refer to this course for a great overview of the Laplacian operator and its use in geometry processing.
For the generalization to volumes, see this pdf
Example
import mouette as M
import numpy as np
mesh = M.mesh.load("path/to/my/mesh/mesh.obj")
W = M.operators.laplacian(mesh,cotan=True)
A = M.operators.area_weight_matrix(mesh, inverse=True)
L = A @ W # discretization of the laplace-beltrami operator is inverted mass matrix times the cotan weight matrix
X = np.random.random(len(mesh.vertices)) # a value in (0,1) per vertex
for _ in range(10):
X = W.dot(X) # perform 10 steps of diffusion
graph_laplacian(mesh)
Simplest Laplacian defined on a graph. uses uniform weights for connectivity.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
mesh
|
Mesh
|
the input mesh |
required |
Returns:
Type | Description |
---|---|
csc_matrix
|
(scipy.sparse.coo_matrix) : the laplacian operator as a sparse matrix |
laplacian(mesh, cotan=True, connection=None, order=4)
Cotan laplacian on vertices.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
mesh
|
SurfaceMesh
|
input mesh |
required |
cotan
|
bool)
|
whether to compute real cotan values for more precise discretization or only 0/1 values as a graph laplacian. Defaults to True. |
True
|
connection
|
SurfaceConnectionVertices
|
For a laplacian on 1-forms, gives the angle in local bases of all adjacent edges. Defaults to None. |
None
|
order
|
int
|
Order of the parallel transport (useful when computing frame fields). Does nothing if parallel_transport is set to None. Defaults to 4. |
4
|
Returns:
Type | Description |
---|---|
csc_matrix
|
scipy.sparse.csc_matrix : the Laplacian operator as a sparse matrix |
laplacian_triangles(mesh, cotan=True, connection=None, order=4)
Cotan laplacian defined on face connectivity (ie on the dual mesh)
Parameters:
Name | Type | Description | Default |
---|---|---|---|
mesh
|
SurfaceMesh
|
the supporting mesh |
required |
cotan
|
bool
|
whether to use cotan laplacian or . Defaults to True. |
True
|
connection
|
SurfaceConnectionFaces
|
description. Defaults to None. |
None
|
order
|
int
|
description. Defaults to 4. |
4
|
Returns:
Type | Description |
---|---|
lil_matrix
|
scipy.sparse.lil_matrix : the Laplacian operator as a sparse matrix |
volume_laplacian(mesh)
Volume laplacian on vertices
This is the 3D extension of the cotan laplacian, ie the discretization of the Laplace-Beltrami operator on 3D manifolds.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
mesh
|
VolumeMesh
|
the input mesh |
required |
Returns:
Type | Description |
---|---|
lil_matrix
|
scipy.sparse.lil_matrix: the Laplacian operator as a sparse matrix |
References
[1] https://cseweb.ucsd.edu/~alchern/projects/ConformalVolume/
[2] https://www.cs.cmu.edu/~kmcrane/Projects/Other/nDCotanFormula.pdf
laplacian_tetrahedra(mesh)
Laplacian defined on cell connectivity (ie on the dual volume mesh)
Parameters:
Name | Type | Description | Default |
---|---|---|---|
mesh
|
SurfaceMesh
|
input mesh |
required |
Returns:
Type | Description |
---|---|
csc_matrix
|
scipy.sparse.csc_matrix |
area_weight_matrix(mesh, inverse=False, sqrt=False, format='csc')
Returns the diagonal matrix A of area weights on vertices.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
mesh
|
SurfaceMesh
|
input mesh |
required |
inverse
|
bool
|
whether to return A or A^-1. Defaults to False. |
False
|
sqrt
|
bool
|
whether to return A or A^{1/2}. Can be combined with |
False
|
format
|
str
|
one of the sparse matrix format of scipy (csc, csr, coo, lil, ...). Defaults to csc. |
'csc'
|
Returns:
Type | Description |
---|---|
csc_matrix
|
sp.csc_matrix: diagonal matrix of vertex areas |
cotan_edge_diagonal(mesh, inverse=True)
Builds a diagonal matrix of size |E|x|E| where the coefficients are the cotan weights, i.e.
M[e,e] = 1/abs( cot(a_e) + cot(b_e))
where a_e and b_e are the opposite angles in adjacent triangles of edge e.
Parameters:
Name | Type | Description | Default |
---|---|---|---|
mesh
|
SurfaceMesh
|
input mesh |
required |
inverse
|
bool
|
whether to compute M or M^-1 (all coefficients on the diagonal inverted) |
True
|
Returns:
Type | Description |
---|---|
csc_matrix
|
sp.csc_matrix |
area_weight_matrix_faces(mesh, inverse=False)
Returns the diagonal matrix A of area weights on faces Laplace-beltrami operator for a 2D manifold is (A^-1)L where A is the area weight and L is the cotan matrix
Parameters:
Name | Type | Description | Default |
---|---|---|---|
mesh
|
SurfaceMesh
|
the input mesh |
required |
inverse
|
bool
|
whether to return A or A^-1. Defaults to False. |
False
|
Returns:
Type | Description |
---|---|
csc_matrix
|
sp.csc_matrix |
volume_weight_matrix(mesh, inverse=False, sqrt=False, format='csc')
Mass diagonal matrix for volume Laplacian.
Args:
mesh (VolumeMesh): input mesh
inverse (bool, optional): whether to return A or A^-1. Defaults to False.
sqrt (bool, optional): whether to return A or A^{1/2}. Can be combined with inverse
to return A^{-1/2}. Defaults to False.
format (str, optional): one of the sparse matrix format of scipy (csc, csr, coo, lil, ...). Defaults to csc.
Returns:
Type | Description |
---|---|
dia_matrix
|
sp.csc_matrix: diagonal matrix of vertices area |