pyFM.mesh.laplacian

Functions

cotangent_weights(vertices, faces)

Compute the cotengenant weights matrix for mesh laplacian.

dia_area_mat(vertices, faces[, faces_areas])

Compute the diagonal matrix of lumped vertex area for mesh laplacian.

fem_area_mat(vertices, faces[, faces_areas])

Compute the area matrix for mesh laplacian using finite elements method.

laplacian_spectrum(W, A[, spectrum_size])

Solves the generalized eigenvalue problem.

pyFM.mesh.laplacian.dia_area_mat(vertices, faces, faces_areas=None)

Compute the diagonal matrix of lumped vertex area for mesh laplacian. Entry i on the diagonal is the area of vertex i, approximated as one third of adjacent triangles

Parameters:
  • vertices – (n,3) array of vertices coordinates

  • faces – (m,3) array of vertex indices defining faces

  • faces_area – (m,) - Optional, array of per-face area

Returns:

A – (n,n) sparse diagonal matrix of vertex areas in dia format

Return type:

scipy.sparse.dia_matrix

pyFM.mesh.laplacian.fem_area_mat(vertices, faces, faces_areas=None)

Compute the area matrix for mesh laplacian using finite elements method.

Entry (i,i) is 1/6 of the sum of the area of surrounding triangles Entry (i,j) is 1/12 of the sum of the area of triangles using edge (i,j)

Parameters:
  • vertices – (n,3) array of vertices coordinates

  • faces – (m,3) array of vertex indices defining faces

  • faces_area – (m,) - Optional, array of per-face area

Returns:

A – (n,n) sparse area matrix in csc format

Return type:

scipy.sparse.csc_matrix

pyFM.mesh.laplacian.cotangent_weights(vertices, faces)

Compute the cotengenant weights matrix for mesh laplacian.

Entry (i,i) is 1/6 of the sum of the area of surrounding triangles Entry (i,j) is 1/12 of the sum of the area of triangles using edge (i,j)

Parameters:
  • vertices – (n,3) array of vertices coordinates

  • faces – (m,3) array of vertex indices defining faces

  • faces_area – (m,) - Optional, array of per-face area

Returns:

A – (n,n) sparse area matrix in csc format

Return type:

scipy.sparse.csc_matrix

pyFM.mesh.laplacian.laplacian_spectrum(W, A, spectrum_size=200)

Solves the generalized eigenvalue problem. Change solver if necessary

Parameters:
  • W – (n,n) - sparse matrix of cotangent weights

  • A – (n,n) - sparse matrix of area weights

  • spectrum_size – int - number of eigenvalues to compute

Returns:

  • eigenvalues (np.ndarray) – (spectrum_size,) - array of eigenvalues

  • eigenvectors (np.ndarray) – (n, spectrum_size) - array of eigenvectors