lapy.Solver¶
- class lapy.Solver(geometry, lump=False, aniso=None, aniso_smooth=10, use_cholmod=False)¶
FEM solver for Laplace Eigenvalue and Poisson Equation.
Inputs can be geometry classes which have vertices and elements. Currently
TriaMesh
andTetMesh
are implemented. FEM matrices (stiffness (or A) and mass matrix (or B)) are computed during the construction. After that the Eigenvalue solver (lapy.Solver.eigs
) or Poisson Solver (lapy.Solver.poisson
) can be called.- Parameters:
- geometry
TriaMesh
|TetMesh
Mesh geometry.
- lump
bool
If True, lump the mass matrix (diagonal).
- aniso
float
|tuple
of shape (2,) Anisotropy for curvature based anisotopic Laplace. If a tuple
(a_min, a_max), differentially affects the minimum and maximum curvature directions. e.g. ``(0, 50)
will set scaling to 1 into the minimum curvature direction, even if the maximum curvature is large in those regions ( i.e. isotropic in regions with large maximum curvature and minimum curvature close to 0, i.e. a concave cylinder).- aniso_smooth
int
|None
Number of smoothing iterations for curvature computation on vertices.
- use_cholmod
bool
, default:False
If True, attempts to use the Cholesky decomposition for improved execution speed. Requires the
scikit-sparse
library. If it can not be found, an error will be thrown. If False, will use slower LU decomposition.
- geometry
Methods
eigs
([k])Compute the linear finite-element method Laplace-Beltrami spectrum.
fem_tria_mass
(tria[, lump])Compute the sparse symmetric mass matrix of the Laplace-Beltrami operator for a given triangle mesh.
poisson
([h, dtup, ntup])Solver for the poisson equation with boundary conditions.
Notes
The class has a static member to create the mass matrix of
TriaMesh
for external function that do not need stiffness.- eigs(k=10)¶
Compute the linear finite-element method Laplace-Beltrami spectrum.
- Parameters:
- k
int
The number of eigenvalues and eigenvectors desired.
k
must be smaller thanN
. It is not possible to compute all eigenvectors of a matrix.
- k
- Returns:
- eigenvalues
array
of shape (k,) Array of k eigenvalues. For closed meshes or Neumann boundary condition,
0
will be the first eigenvalue (with constant eigenvector).- eigenvectors
array
of shape (N, k) Array representing the k eigenvectors. The column
eigenvectors[:, i]
is the eigenvector corresponding toeigenvalues[i]
.
- eigenvalues
- static fem_tria_mass(tria, lump=False)¶
Compute the sparse symmetric mass matrix of the Laplace-Beltrami operator for a given triangle mesh.
The sparse symmetric matrix is computed for a given triangle mesh using the linear finite element method (assuming a closed mesh or Neumann boundary condition). This function is faster than the constructor above and can be used when only a mass matrix is needed.
- Parameters:
- Returns:
- B
csc_matrix
of shape (n, n) Sparse symmetric positive definite matrix.
- B
Notes
This only returns the mass matrix
B
of the Eigenvalue problem:A x = lambda B x
. The area of the surface mesh can be obtained viaB.sum()
.
- poisson(h=0.0, dtup=(), ntup=())¶
Solver for the poisson equation with boundary conditions.
This solver is based on the
A
andB
Laplace matrices whereA x = B h
andA
is a sparse symmetric positive semi definitive matrix of shape(n`, n)
and B is a sparse symmetric positive definitive matrix of shape(n, n)
.- Parameters:
- h
float
|array
Right hand side, can be constant or array with vertex values. The default
0
corresponds to Laplace equationA x = 0
.- dtup
tuple
Dirichlet boundary condition as a tuple containing the index and data arrays of same length. The default, an empty tuple, corresponds to no Dirichlet condition.
- ntup
tuple
Neumann boundary condition as a tuple containing the index and data arrays of same length. The default, an empty tuple, corresponds to Neumann on all boundaries.
- h
- Returns:
- x
array
Array with vertex value of the solution.
- x
Notes
A
andB
are obtained viacomputeAB
for either triangle or tetraheral mesh.