pyFM.mesh.trimesh

Classes

TriMesh(*args, **kwargs)

Mesh (and PointCloud) Class

class pyFM.mesh.trimesh.TriMesh(*args, **kwargs)

Mesh (and PointCloud) Class

Parameters:
  • path (str, optional) – path to a .off file

  • vertices (np.ndarray, optional) – (n,3) vertices coordinates

  • faces (np.ndarray, optional) – (m,3) list of indices of triangles

  • area_normalize (bool, optional) – If True, normalize the mesh

  • center (bool, optional) – If True, center the mesh

  • rotation (np.ndarray, optional) – 3x3 rotation matrix

  • translation (np.ndarray, optional) – 3D translation vector, applied after rotation

path

path the the loaded .off file. Set to None if the geometry is modified.

Type:

str

meshname
name of the .off file. Remains even when geometry is modified. ‘_n’ is

added at the end if the mesh was normalized.

Type:

str

W

(n,n) sparse cotangent weight matrix

A

(n,n) sparse area matrix (either diagonal or computed with finite elements)

eigenvalues

(K,) eigenvalues of the Laplace Beltrami Operator

eigenvectors

(n,K) eigenvectors of the Laplace Beltrami Operator

property vertlist

Get or set the vertices. Checks the format when setting

Returns:

vertlist – (n,3) array of vertices

Return type:

np.ndarray

property facelist

Get or set the faces. Checks the format when setting

Returns:

facelist – (m,3) array of faces

Return type:

np.ndarray

property vertices

alias for vertlist

Returns:

vertices – (n,3) array of vertices

Return type:

np.ndarray

property faces

alias for facelist

Returns:

faces – (m,3) array of faces

Return type:

np.ndarray

property n_vertices

return the number of vertices in the mesh

Returns:

n_vertices – number of vertices in the mesh

Return type:

int

property n_faces

return the number of faces in the mesh

Returns:

n_faces – number of faces in the mesh

Return type:

int

property area

Returns the area of the mesh

Returns:

area – area of the mesh

Return type:

float

property sqrtarea

square root of the area

Returns:

sqrtarea – square root of the area

Return type:

float

property edges

return a (p,2) array of edges defined by vertex indices.

Returns:

edges – (p,2) array of edges

Return type:

np.ndarray

property normals

return face normals

Returns:

normals – (m,3) array of face normals

Return type:

np.ndarray

property vertex_normals

Returns per vertex_normal

Returns:

vertex_normals – (n,3) array of vertex normals

Return type:

np.ndarray

property vertex_areas

per vertex area

Returns:

vertex_areas – (n,) array of vertex areas

Return type:

np.ndarray

property faces_areas

per face area

Returns:

faces_areas – (m,) array of face areas

Return type:

np.ndarray

property face_areas

per face area

Returns:

faces_areas – (m,) array of face areas

Return type:

np.ndarray

property center_mass

center of mass

Returns:

center_mass – (3,) array of the center of mass

Return type:

np.ndarray

property is_normalized

Whether the mash has been manually normalized using the self.area_normalize method

Returns:

is_normalized – Whether the mesh has been area normalized

Return type:

bool

property is_modified

Whether the mash has been modified from path with non-isometric deformations

Returns:

is_modified – Whether the mesh has been modified wrt to original input

Return type:

bool

area_normalize()

Normalize the mesh by its area

rotate(R)

Rotate mesh and normals

Parameters:

R (np.ndarray) – (3,3) rotation matrix

translate(t)

translate mesh

Parameters:

t (np.ndarray) – (3,) translation vector

scale(alpha)

Multiply mesh by alpha. modify vertices, area, spectrum, geodesic distances

Parameters:

alpha (float) – scaling factor

center()

center the mesh

laplacian_spectrum(k, intrinsic=False, return_spectrum=True, robust=False, verbose=False)

Compute the Laplace Beltrami Operator and its spectrum. Consider using the .process() function for easier use !

Parameters:
  • K (int) – number of eigenvalues to compute

  • intrinsic (bool, optional) – Use intrinsic triangulation. Defaults to false

  • robust (bool, optional) – use tufted laplacian, defaults to False

  • return_spectrum (bool, optional) – Whether to return the computed spectrum, defaults to True

Returns:

  • eigenvalues (np.ndarray, optional) – (k,) - Only if return_spectrum is True.

  • eigenvectors (np.ndarray, optional) – (n,k) - Only if return_spectrum is True.

process(k=200, skip_normals=True, intrinsic=False, robust=False, verbose=False)

Process the LB spectrum and saves it. Additionnaly computes per-face normals

Parameters:
  • k (int) – (default = 200) Number of eigenvalues to compute

  • skip_normals (bool, optional) – If set to True, skip normals computation. Defaults to True

  • intrinsic (bool, optional) – Use intrinsic triangulation. Defaults to False

  • robust (bool) – use tufted laplacian

  • verbose (bool) – print progress

project(func, k=None)

Project one or multiple functions on the spectral basis

Parameters:
  • func (np.ndarray) – (n,p) or (n,) functions on the shape

  • k (int) – dimension of the LB basis on which to project. If None use all the computed basis

Returns:

projected_func – (k,p) or (k,) projected function

Return type:

np.ndarray

decode(projection)

Build a function from its coefficient in the spectral basis

Parameters:

projection (np.ndarray) – (k,p) or (k,) functions on the reduced basis of the shape

Returns:

func – (n,p) or (n,) projected function

Return type:

np.ndarray

unproject(projection)

Alias for decode

Parameters:

projection (np.ndarray) – (k,p) or (k,) functions on the reduced basis of the shape

reconstruct(func, k=None)

Reconstruct function with the LB eigenbasis, ie project on the spectral basis and rebuild values on all vertices.

Parameters:
  • func (np.ndarray) – (n,p) or (n,) - functions on the shape

  • k (int) – Number of eigenfunctions to use. If None, uses the complete computed basis.

Returns:

func – (n,p) or (n,) projected function

Return type:

np.ndarray

get_geodesic(dijkstra=False, robust=True, save=False, force_compute=False, sym=False, batch_size=500, verbose=False)

Compute the geodesic distance matrix using either the Dijkstra algorithm or the Heat Method. Loads from cache if possible.

Parameters:
  • dijkstra (bool , optional) – If True, use Dijkstra algorithm instead of the heat method. Defaults to False

  • robust (boo, optional) – Robust heat method. Defaults to True

  • save (bool, optional) – If True, save the resulting distance matrix at ‘{path}/geod_cache/{meshname}.npy’ with ‘path/meshname.{ext}’ path of the current mesh. Defaults to False

  • force_compute (bool, optional) – If True, doesn’t look for a cached distance matrix. Defaults to False

  • sym (bool, optional) – Symmetrize the matrix if computed with heat method. Defaults to False

  • batch_size (int, optional) – If robust is False, compute distances by batch

  • verbose (bool, optional) – Print progress

Returns:

distances – (n,n) matrix of geodesic distances

Return type:

np.ndarray

geod_from(i, robust=True)

Compute geodesic distances from vertex i sing the Heat Method

Parameters:
  • i (int) – index from source

  • robust (bool, optional) – Robust heat method

Returns:

dist – (n,) distances to vertex i

Return type:

np.ndarray

l2_sqnorm(func)

Return the squared L2 norm of one or multiple functions on the mesh. For a single function f, this returns f.T @ A @ f with A the area matrix.

Parameters:

func (np.ndarray) – (n,p) or (n,) functions on the mesh

Returns:

sqnorm – (p,) array of squared l2 norms or a float only one function was provided.

Return type:

np.ndarray

l2_inner(func1, func2)

Return the L2 inner product of two functions, or pairwise inner products if lists of function is given.

For two functions f1 and f2, this returns f1.T @ A @ f2 with A the area matrix.

Parameters:
  • func1 (np.ndarray) – (n,p) or (n,) functions on the mesh

  • func2 (np.ndarray) – (n,p) or (n,) functions on the mesh

Returns:

sqnorm

(p,) array of L2 inner product or a float only one function per argument

was provided.

Return type:

np.ndarray

h1_sqnorm(func)

Return the squared H^1_0 norm (L2 norm of the gradient) of one or multiple functions on the mesh. For a single function f, this returns f.T @ W @ f with W the stiffness matrix.

Parameters:

func (np.ndarray) – (n,p) or (n,) functions on the mesh

Returns:

sqnorm – (p,) array of squared H1 norms or a float only one function was provided.

Return type:

np.ndarray

h1_inner(func1, func2)

Return the H1 inner product of two functions, or pairwise inner products if lists of function is given.

For two functions f1 and f2, this returns f1.T @ W @ f2 with W the stiffness matrix.

Parameters:
  • func1 (np.ndarray) – (n,p) or (n,) functions on the mesh

  • func2 (np.ndarray) – (n,p) or (n,) functions on the mesh

Returns:

sqnorm

(p,) array of H1 inner product or a float only one function per argument

was provided.

Return type:

np.ndarray

integrate(func)

Integrate a function or a set of function on the mesh

Parameters:

func (np.ndarray) – (n,p) or (n,) functions on the mesh

Returns:

integral – (p,) array of integrals or a float only one function was provided.

Return type:

np.ndarray

extract_fps(size, random_init=True, geodesic=True, no_load=False, verbose=False)

Samples points using farthest point sampling with geodesic distances. If the geodesic matrix is precomputed (in the cache folder) uses it, else computes geodesic distance in real time

Parameters:
  • size (int) – number of points to sample

  • random_init (bool, optional) – Whether to sample the first point randomly or to take the furthest away from all the other ones. This is only done if the geodesic matrix is accessible from cache. defaults to True

  • geodesic (bool, optional) – If True perform geodesic fps, else euclidean. Defaults to True

  • no_load (bool, optional) – if True never loads cache. Defaults to False

  • verbose (bool, optional) – Print progress. Defaults to False

Returns:

fps – (size,) array of indices of sampled points (given on the complete mesh)

Return type:

np.ndarray

extract_fps_sub(size, sub_points, return_sub_inds=False, random_init=True, geodesic=True, no_load=False, verbose=False)

Samples points using farthest point sampling with geodesic distances, but reduced on a set of samples. If the geodesic matrix is precomputed (in the cache folder) uses it, else computes geodesic distance in real time

Parameters:
  • size (int) – number of points to sample

  • sub_points (np.ndarray) – (size,) array of indices of the sub points

  • random_init – Whether to sample the first point randomly or to take the furthest away from all the other ones. This is only done if the geodesic matrix is accessible from cache. defaults to True

  • geodesic (bool, optional) – If True perform geodesic fps, else eucliden. Defaults to True

  • no_load (bool) – if True never loads cache. Defaults to False

  • verbose (bool) – Print progress. Defaults to False

Returns:

  • fps (np.ndarray) – (size,) array of indices of sampled points (given on the complete mesh)

  • fps_sub (np.ndarray) – (size,) array of indices of sampled points (given on the sub mesh)

gradient(f, normalize=False)

computes the gradient of a function on f using linear interpolation between vertices.

Parameters:
  • f (np.ndarray) – (n_v,) function value on each vertex

  • normalize (bool) – Whether the gradient should be normalized on each face

Returns:

gradient – (n_f,3) gradient of f on each face

Return type:

np.ndarray

divergence(f)

Computes the divergence of a vector field on the mesh

Parameters:

f (np.ndarray) – (n_f, 3) vector value on each face

Returns:

divergence – (n_v,) divergence of f on each vertex

Return type:

np.ndarray

orientation_op(gradf, normalize=False)

Compute the orientation operator associated to a gradient field gradf.

For a given function g on the vertices, this operator linearly computes < grad(f) x grad(g), n> for each vertex by averaging along the adjacent faces. In practice, we compute < n x grad(f), grad(g) > for simpler computation.

Parameters:
  • gradf (np.ndarray) – (n_f,3) gradient field on the mesh

  • normalize (bool, optional) – Whether to normalize the gradient on each face

Returns:

operator – (n_v,n_v) orientation operator.

Return type:

scipy.sparse.csc_matrix

export(filename, precision=None)

Write the mesh in a .off file

Parameters:
  • filename (str) – path to the file to write

  • precision (int) – floating point precision

get_uv(ind1, ind2, mult_const, rotation=None)

Extracts UV coordinates for each vertices.

Extracted by orthogonal projection on 2 of x,y,z axes

Parameters:
  • ind1 (int) – column index to use as first coordinate

  • ind2 (int) – column index to use as second coordinate

  • mult_const (float) – number of time to repeat the pattern

Returns:

uv – (n,2) UV coordinates of each vertex

Return type:

np.ndarray

export_texture(filename, uv, mtl_file='material.mtl', texture_im='texture_1.jpg', precision=None, verbose=False)

Write a .obj file with texture using uv coordinates

Parameters:
  • filename (str) – path to the .obj file to write

  • uv – (n,2) uv coordinates of each vertex

  • mtl_file (str) – name of the .mtl file

  • texture_im (str) – name of the .jpg file definig texture

compute_normals()

Compute normal vectors for each face

Returns:

normals – (m,3) array of normal vectors

Return type:

np.ndarray

set_vertex_normal_weighting(weight_type)

Set weighting type for vertex normals between ‘area’ and ‘uniform’

compute_vertex_normals()

computes vertex normals in self.vertex_normals

Returns:

vertex_normals – (n,3) array of vertex normals

Return type:

np.ndarray

compute_edges()

computes edges in self.edges

Returns:

edges – (p,2) array of edges

Return type:

np.ndarray