quimb.linalg.slepc_linalg¶
Interface to slepc4py for solving advanced eigenvalue problems.
Attributes¶
Classes¶
Helper class for caching the PETSc and SLEPc imports and making sure |
|
Functions¶
Define the default communicator. |
|
|
|
|
|
|
|
|
Slice the matrix mat between indices ri and rf -- for csr or bsr. |
|
Convert a matrix to the relevant PETSc type, currently |
|
Convert a vector/array to the PETSc form. |
|
Create an empty petsc vector of size d. |
|
Gather the petsc vector/matrix x to a single array on the master |
|
Take the real part of a set of vectors and normalize it. This is used |
|
Initialise a krylov subspace and preconditioner. |
|
Create a slepc spectral transformation object with specified solver. |
|
|
|
Create an advanced eigensystem solver |
|
Solve a matrix using the advanced eigensystem solver |
|
|
|
Find the singular values for sparse matrix a. |
|
Compute the action of |
Look up PETSc error to print when raising after not converging. |
|
|
Module Contents¶
- class quimb.linalg.slepc_linalg._PetscSlepcHandler[source]¶
Helper class for caching the PETSc and SLEPc imports and making sure that the MPI comm they use is aligned.
- _comm = '__UNINITIALIZED__'¶
- quimb.linalg.slepc_linalg._petsc_slepc_handler¶
- quimb.linalg.slepc_linalg.slice_sparse_matrix_to_components(mat, ri, rf)[source]¶
Slice the matrix mat between indices ri and rf – for csr or bsr.
- quimb.linalg.slepc_linalg.convert_mat_to_petsc(mat, comm=None)[source]¶
Convert a matrix to the relevant PETSc type, currently only supports csr, bsr and dense matrices formats.
- Parameters:
mat (dense, sparse, LinearOperator or Lazy matrix.) – The operator to convert.
comm (mpi4py.MPI.Comm instance) – The mpi communicator.
- Returns:
pmat – The matrix in petsc form - only the local part if running across several mpi processes.
- Return type:
petsc4py.PETSc.Mat
- quimb.linalg.slepc_linalg.convert_vec_to_petsc(vec, comm=None)[source]¶
Convert a vector/array to the PETSc form.
- Parameters:
vec (vector-like) – Numpy array, will be unravelled to one dimension.
comm (mpi4py.MPI.Comm instance) – The mpi communicator.
- Returns:
pvec – The vector in petsc form - only the local part if running across several mpi processes.
- Return type:
petsc4py.PETSc.Vec
- quimb.linalg.slepc_linalg.new_petsc_vec(d, comm=None)[source]¶
Create an empty petsc vector of size d.
- Parameters:
d (int) – Dimension of vector, i.e. the global size.
comm (mpi4py.MPI.Comm instance) – The mpi communicator.
- Returns:
pvec – An empty vector in petsc form - only the local part if running across several mpi processes.
- Return type:
petsc4py.PETSc.Vec
- quimb.linalg.slepc_linalg.gather_petsc_array(x, comm, out_shape=None)[source]¶
Gather the petsc vector/matrix x to a single array on the master process, assuming that owernership is sliced along the first dimension.
- Parameters:
x (petsc4py.PETSc Mat or Vec) – Distributed petsc array to gather.
comm (mpi4py.MPI.COMM) – MPI communicator
out_shape (tuple, optional) – If not None, reshape the output array to this.
- Returns:
gathered
- Return type:
np.array master, None on workers (rank > 0)
- quimb.linalg.slepc_linalg.normalize_real_part(vecs, transposed=False)[source]¶
Take the real part of a set of vectors and normalize it. This is used for returning real eigenvectors even when the PETSc scalar type is complex.
- quimb.linalg.slepc_linalg._init_krylov_subspace(comm=None, tol=None, maxiter=None, KSPType='preonly', PCType='lu', PCFactorSolverType='mumps')[source]¶
Initialise a krylov subspace and preconditioner.
- quimb.linalg.slepc_linalg._init_spectral_inverter(STType='sinvert', KSPType='preonly', PCType='lu', PCFactorSolverType='mumps', comm=None)[source]¶
Create a slepc spectral transformation object with specified solver.
- quimb.linalg.slepc_linalg._WHICH_SCIPY_TO_SLEPC¶
- quimb.linalg.slepc_linalg._init_eigensolver(k=6, which='LM', sigma=None, isherm=True, isgen=False, EPSType=None, st_opts=None, tol=None, A_is_linop=False, B_is_linop=False, maxiter=None, ncv=None, l_win=None, comm=None)[source]¶
Create an advanced eigensystem solver
- Parameters:
sigma – Target eigenvalue.
isherm – Whether problem is hermitian or not.
- Return type:
SLEPc solver ready to be called.
- quimb.linalg.slepc_linalg.eigs_slepc(A, k, *, B=None, which=None, sigma=None, isherm=True, P=None, v0=None, ncv=None, return_vecs=True, sort=True, EPSType=None, return_all_conv=False, st_opts=None, tol=None, maxiter=None, l_win=None, comm=None)[source]¶
Solve a matrix using the advanced eigensystem solver
- Parameters:
A (dense-matrix, sparse-matrix, LinearOperator or callable) – Operator to solve.
k (int, optional) – Number of requested eigenpairs.
B (dense-matrix, sparse-matrix, LinearOperator or callable) – The RHS operator defining a generalized eigenproblem.
which ({"LM": "SM", "LR", "LA", "SR", "SA", "LI", "SI", "TM", "TR", "TI"}) – Which eigenpairs to target. See
scipy.sparse.linalg.eigs()
.sigma (float, optional) – Target eigenvalue, implies
which='TR'
if this is not set butsigma
is.isherm (bool, optional) – Whether problem is hermitian or not.
P (dense-matrix, sparse-matrix, LinearOperator or callable) – Perform the eigensolve in the subspace defined by this projector.
v0 (1D-array like, optional) – Initial iteration vector, e.g., informed guess at eigenvector.
ncv (int, optional) – Subspace size, defaults to
min(20, 2 * k)
.return_vecs (bool, optional) – Whether to return the eigenvectors.
sort (bool, optional) – Whether to sort the eigenpairs in ascending real value.
EPSType ({"krylovschur", 'gd', 'lobpcg', 'jd', 'ciss'}, optional) – SLEPc eigensolver type to use, see slepc4py.EPSType.
return_all_conv (bool, optional) – Whether to return converged eigenpairs beyond requested subspace size
st_opts (dict, optional) – options to send to the eigensolver internal inverter.
tol (float, optional) – Tolerance.
maxiter (int, optional) – Maximum number of iterations.
comm (mpi4py communicator, optional) – MPI communicator, defaults to
COMM_SELF
for a single process solve.
- Returns:
lk ((k,) array) – The eigenvalues.
vk ((m, k) array) – Corresponding eigenvectors (if
return_vecs=True
).
- quimb.linalg.slepc_linalg._init_svd_solver(nsv=6, SVDType='cross', tol=None, maxiter=None, ncv=None, comm=None)[source]¶
- quimb.linalg.slepc_linalg.svds_slepc(A, k=6, ncv=None, return_vecs=True, SVDType='cross', return_all_conv=False, tol=None, maxiter=None, comm=None)[source]¶
Find the singular values for sparse matrix a.
- Parameters:
A (sparse matrix in csr format) – The matrix to solve.
k (int) – Number of requested singular values.
method ({"cross", "cyclic", "lanczos", "trlanczos"}) – Solver method to use.
- Returns:
U ((m, k) array) – Left singular vectors (if
return_vecs=True
) as columns.s ((k,) array) – Singular values.
VH ((k, n) array) – Right singular vectors (if
return_vecs=True
) as rows.
- quimb.linalg.slepc_linalg.mfn_multiply_slepc(mat, vec, fntype='exp', MFNType='AUTO', comm=None, isherm=False)[source]¶
Compute the action of
func(mat) @ vec
.- Parameters:
mat (operator) – Operator to compute function action of.
vec (vector-like) – Vector to compute matrix function action on.
func ({'exp', 'sqrt', 'log'}, optional) – Function to use.
MFNType ({'krylov', 'expokit'}, optional) – Method of computing the matrix function action, ‘expokit’ is only available for func=’exp’.
comm (mpi4py.MPI.Comm instance, optional) – The mpi communicator.
isherm (bool, optional) – If mat is known to be hermitian, this might speed things up in some circumstances.
- Returns:
fvec – The vector output of
func(mat) @ vec
.- Return type:
array