slepc4py.SLEPc.EPS#
- class slepc4py.SLEPc.EPS#
Bases:
ObjectEigenvalue Problem Solver.
The Eigenvalue Problem Solver (
EPS) is the object provided by slepc4py for specifying a linear eigenvalue problem, either in standard or generalized form. It provides uniform and efficient access to all of the linear eigensolvers included in the package.Enumerations
EPS type of balancing used for non-Hermitian problems.
EPS CISS extraction technique.
EPS CISS quadrature rule.
EPS convergence test.
EPS convergence reasons.
EPS error type to assess accuracy of computed solutions.
EPS extraction technique.
EPS Krylov-Schur method for BSE problems.
EPS Lanczos reorthogonalization type.
EPS Power shift type.
EPS problem type.
EPS stopping test.
EPS type.
EPS desired part of spectrum.
Methods Summary
appendOptionsPrefix([prefix])Append to the prefix used for searching for all EPS options in the database.
Clear all monitors for an
EPSobject.computeError(i[, etype])Compute the error associated with the i-th computed eigenpair.
create([comm])Create the EPS object.
destroy()Destroy the EPS object.
errorView([etype, viewer])Display the errors associated with the computed solution.
Get the arbitrary selection function.
Get the type of reorthogonalization used during the Arnoldi iteration.
getBV()Get the basis vectors object associated to the eigensolver.
Get the balancing type used by the EPS, and the associated parameters.
Get the extraction technique used in the CISS solver.
Get the array of linear solver objects associated with the CISS solver.
Get the quadrature rule used in the CISS solver.
Get the values of various refinement parameters in the CISS solver.
Get the values of various size parameters in the CISS solver.
Get the values of various threshold parameters in the CISS solver.
Get the flag indicating the use of the
STobject in the CISS solver.Get the number of converged eigenpairs.
Get the reason why the
solve()iteration was stopped.Get how to compute the error estimate used in the convergence test.
getDS()Get the direct solver associated to the eigensolver.
Get number of eigenvalues to compute and the dimension of the subspace.
getEigenpair(i[, Vr, Vi])Get the i-th solution of the eigenproblem as computed by
solve().Get the i-th eigenvalue as computed by
solve().Get the eigenvalue comparison function.
getEigenvector(i[, Vr, Vi])Get the i-th right eigenvector as computed by
solve().Get the error estimate associated to the i-th computed eigenpair.
Get the extraction type used by the EPS object.
Get the orthogonalization used in the search subspace.
Get the number of vectors to be added to the searching space.
Get a flag indicating whether the double expansion variant is active.
Get the initial size of the searching space.
Get a flag indicating if the search subspace is started with a Krylov basis.
Get the number of vectors of the search space after restart.
Get the computational interval for spectrum slicing.
Get an orthonormal basis of the computed invariant subspace.
Get the current iteration number.
Get the orthogonalization used in the search subspace.
Get the number of vectors to be added to the searching space.
Get the flag indicating if the dynamic stopping is being used.
getJDFix()Get the threshold for changing the target in the correction equation.
Get the initial size of the searching space.
Get a flag indicating if the search subspace is started with a Krylov basis.
Get the number of vectors of the search space after restart.
Get the method used for BSE structured eigenproblems (Krylov-Schur).
Get the flag that enforces zero detection in spectrum slicing.
Get the dimensions used for each subsolve step (spectrum slicing).
Get the values of the shifts and their corresponding inertias.
Get the linear solver object associated with the internal
EPSobject.Get the locking flag used in the Krylov-Schur method.
Get the number of partitions of the communicator (spectrum slicing).
Get the restart parameter used in the Krylov-Schur method.
Get information related to the case of doing spectrum slicing.
Get the eigenproblem matrices stored in the subcommunicator.
getKrylovSchurSubcommPairs(i[, v])Get the i-th eigenpair stored in the multi-communicator of the process.
Get the points that delimit the subintervals.
Get the block size used in the LOBPCG method.
Get the locking flag used in the LOBPCG method.
Get the restart parameter used in the LOBPCG method.
Get the type of reorthogonalization used during the Lanczos iteration.
getLeftEigenvector(i[, Wr, Wi])Get the i-th left eigenvector as computed by
solve().Get the rank values used for the Lyapunov step.
Get the list of monitor functions.
Get the matrices associated with the eigenvalue problem.
Get the prefix used for searching for all EPS options in the database.
Get the type of shifts used during the power iteration.
Get the problem type from the EPS object.
Get the flag indicating whether purification is activated or not.
getRG()Get the region object associated to the eigensolver.
Get the reset parameter used in the RQCG method.
getST()Get the spectral transformation object associated to the eigensolver.
Get the stopping test function.
Get the value of the target.
Get the threshold used in the threshold stopping test.
Get the tolerance and max.
Get the flag indicating if all residual norms must be computed or not.
Get the flag indicating if true residual must be computed explicitly.
Get the flag indicating if a two-sided variant of the algorithm is being used.
getType()Get the EPS type of this object.
Get which portion of the spectrum is to be sought.
Tell if the EPS object corresponds to a generalized eigenproblem.
Tell if the EPS object corresponds to a Hermitian eigenproblem.
Eigenproblem requiring a positive (semi-) definite matrix \(B\).
Tell if the EPS object corresponds to a structured eigenvalue problem.
reset()Reset the EPS object.
setArbitrarySelection(arbitrary[, args, kargs])Set an arbitrary selection criterion function.
setArnoldiDelayed(delayed)Set (toggle) delayed reorthogonalization in the Arnoldi iteration.
setBV(bv)Set a basis vectors object associated to the eigensolver.
setBalance([balance, iterations, cutoff])Set the balancing technique to be used by the eigensolver.
setCISSExtraction(extraction)Set the extraction technique used in the CISS solver.
setCISSQuadRule(quad)Set the quadrature rule used in the CISS solver.
setCISSRefinement([inner, blsize])Set the values of various refinement parameters in the CISS solver.
setCISSSizes([ip, bs, ms, npart, bsmax, ...])Set the values of various size parameters in the CISS solver.
setCISSThreshold([delta, spur])Set the values of various threshold parameters in the CISS solver.
setCISSUseST(usest)Set a flag indicating that the CISS solver will use the
STobject.setConvergenceTest(conv)Set how to compute the error estimate used in the convergence test.
setDS(ds)Set a direct solver object associated to the eigensolver.
setDeflationSpace(space)Set vectors to form a basis of the deflation space.
setDimensions([nev, ncv, mpd])Set number of eigenvalues to compute and the dimension of the subspace.
setEigenvalueComparison(comparison[, args, ...])Set an eigenvalue comparison function.
setExtraction(extraction)Set the extraction type used by the eigensolver.
Set EPS options from the options database.
setGDBOrth(borth)Set the orthogonalization that will be used in the search subspace.
setGDBlockSize(bs)Set the number of vectors to be added to the searching space.
setGDDoubleExpansion(doubleexp)Set that the search subspace is expanded with double expansion.
setGDInitialSize(initialsize)Set the initial size of the searching space.
setGDKrylovStart([krylovstart])Set (toggle) starting the search subspace with a Krylov basis.
setGDRestart([minv, plusk])Set the number of vectors of the search space after restart.
setInitialSpace(space)Set the initial space from which the eigensolver starts to iterate.
setInterval(inta, intb)Set the computational interval for spectrum slicing.
setJDBOrth(borth)Set the orthogonalization that will be used in the search subspace.
setJDBlockSize(bs)Set the number of vectors to be added to the searching space.
setJDConstCorrectionTol(constant)Deactivate the dynamic stopping criterion.
setJDFix(fix)Set the threshold for changing the target in the correction equation.
setJDInitialSize(initialsize)Set the initial size of the searching space.
setJDKrylovStart([krylovstart])Set (toggle) starting the search subspace with a Krylov basis.
setJDRestart([minv, plusk])Set the number of vectors of the search space after restart.
Set the Krylov-Schur variant used for BSE structured eigenproblems.
setKrylovSchurDetectZeros(detect)Set the flag that enforces zero detection in spectrum slicing.
setKrylovSchurDimensions([nev, ncv, mpd])Set the dimensions used for each subsolve step (spectrum slicing).
setKrylovSchurLocking(lock)Set (toggle) locking/non-locking variants of the Krylov-Schur method.
setKrylovSchurPartitions(npart)Set the number of partitions of the communicator (spectrum slicing).
setKrylovSchurRestart(keep)Set the restart parameter for the Krylov-Schur method.
setKrylovSchurSubintervals(subint)Set the subinterval boundaries.
Set the block size of the LOBPCG method.
setLOBPCGLocking(lock)Toggle between locking and non-locking (LOBPCG method).
setLOBPCGRestart(restart)Set the restart parameter for the LOBPCG method.
setLanczosReorthogType(reorthog)Set the type of reorthogonalization used during the Lanczos iteration.
setLeftInitialSpace(space)Set a left initial space from which the eigensolver starts to iterate.
setLyapIIRanks([rkc, rkl])Set the ranks used in the solution of the Lyapunov equation.
setMonitor(monitor[, args, kargs])Append a monitor function to the list of monitors.
setOperators(A[, B])Set the matrices associated with the eigenvalue problem.
setOptionsPrefix([prefix])Set the prefix used for searching for all EPS options in the database.
setPowerShiftType(shift)Set the type of shifts used during the power iteration.
setProblemType(problem_type)Set the type of the eigenvalue problem.
setPurify([purify])Set (toggle) eigenvector purification.
setRG(rg)Set a region object associated to the eigensolver.
setRQCGReset(nrest)Set the reset parameter of the RQCG iteration.
setST(st)Set a spectral transformation object associated to the eigensolver.
setStoppingTest(stopping[, args, kargs])Set a function to decide when to stop the outer iteration of the eigensolver.
setTarget(target)Set the value of the target.
setThreshold(thres[, rel])Set the threshold used in the threshold stopping test.
setTolerances([tol, max_it])Set the tolerance and max.
setTrackAll(trackall)Set if the solver must compute the residual of all approximate eigenpairs.
setTrueResidual(trueres)Set if the solver must compute the true residual explicitly or not.
setTwoSided(twosided)Set to use a two-sided variant that also computes left eigenvectors.
setType(eps_type)Set the particular solver to be used in the EPS object.
setUp()Set up all the internal data structures.
setWhichEigenpairs(which)Set which portion of the spectrum is to be sought.
solve()Solve the eigensystem.
updateKrylovSchurSubcommMats([s, a, Au, t, ...])Update the eigenproblem matrices stored internally in the communicator.
valuesView([viewer])Display the computed eigenvalues in a viewer.
vectorsView([viewer])Output computed eigenvectors to a viewer.
view([viewer])Print the EPS data structure.
Attributes Summary
The basis vectors (
BV) object associated.The direct solver (
DS) object associated.The type of extraction technique to be employed.
The maximum iteration count.
The type of the eigenvalue problem.
Eigenvector purification.
The region (
RG) object associated.The spectral transformation (
ST) object associated.The value of the target.
The tolerance.
Compute the residual norm of all approximate eigenpairs.
Compute the true residual explicitly.
Two-sided that also computes left eigenvectors.
The portion of the spectrum to be sought.
Methods Documentation
- appendOptionsPrefix(prefix=None)#
Append to the prefix used for searching for all EPS options in the database.
Logically collective.
- cancelMonitor()#
Clear all monitors for an
EPSobject.Logically collective.
See also
Source code at slepc4py/SLEPc/EPS.pyx:1860
- Return type:
- computeError(i, etype=None)#
Compute the error associated with the i-th computed eigenpair.
Collective.
Compute the error (based on the residual norm) associated with the i-th computed eigenpair.
- Parameters:
- Returns:
The error bound, computed in various ways from the residual norm \(\|Ax-\lambda Bx\|_2\) where \(\lambda\) is the eigenvalue and \(x\) is the eigenvector.
- Return type:
Notes
The index
ishould be a value between0andnconv-1(seegetConverged()).If the computation of left eigenvectors was enabled with
setTwoSided(), then the error will be computed using the maximum of the value above and the left residual norm \(\|y^*A-\lambda y^*B\|_2\), where \(y\) is the approximate left eigenvector.See also
- create(comm=None)#
Create the EPS object.
Collective.
- Parameters:
comm (Comm | None) – MPI communicator; if not provided, it defaults to all processes.
- Return type:
See also
- destroy()#
Destroy the EPS object.
Collective.
See also
Source code at slepc4py/SLEPc/EPS.pyx:370
- Return type:
- errorView(etype=None, viewer=None)#
Display the errors associated with the computed solution.
Collective.
Display the errors and the eigenvalues.
- Parameters:
viewer (petsc4py.PETSc.Viewer | None) – Visualization context; if not provided, the standard output is used.
- Return type:
Notes
By default, this function checks the error of all eigenpairs and prints the eigenvalues if all of them are below the requested tolerance. If the viewer has format
ASCII_INFO_DETAILthen a table with eigenvalues and corresponding errors is printed.See also
- getArbitrarySelection()#
Get the arbitrary selection function.
Not collective.
- Returns:
The arbitrary selection function.
- Return type:
See also
- getArnoldiDelayed()#
Get the type of reorthogonalization used during the Arnoldi iteration.
Not collective.
- Returns:
Trueif delayed reorthogonalization is to be used.- Return type:
See also
- getBV()#
Get the basis vectors object associated to the eigensolver.
Not collective.
- Returns:
The basis vectors context.
- Return type:
- getBalance()#
Get the balancing type used by the EPS, and the associated parameters.
Not collective.
- Returns:
- Return type:
See also
- getCISSExtraction()#
Get the extraction technique used in the CISS solver.
Not collective.
- Returns:
The extraction technique.
- Return type:
See also
- getCISSKSPs()#
Get the array of linear solver objects associated with the CISS solver.
Not collective.
- Returns:
The linear solver objects.
- Return type:
Notes
The number of
petsc4py.PETSc.KSPsolvers is equal to the number of integration points divided by the number of partitions. This value is halved in the case of real matrices with a region centered at the real axis.See also
- getCISSQuadRule()#
Get the quadrature rule used in the CISS solver.
Not collective.
- Returns:
The quadrature rule.
- Return type:
See also
- getCISSRefinement()#
Get the values of various refinement parameters in the CISS solver.
Not collective.
- Returns:
- Return type:
See also
- getCISSSizes()#
Get the values of various size parameters in the CISS solver.
Not collective.
- Returns:
- Return type:
See also
- getCISSThreshold()#
Get the values of various threshold parameters in the CISS solver.
Not collective.
- Returns:
- Return type:
See also
- getCISSUseST()#
Get the flag indicating the use of the
STobject in the CISS solver.Not collective.
See also
- getConverged()#
Get the number of converged eigenpairs.
Not collective.
- Returns:
nconv – Number of converged eigenpairs.
- Return type:
Notes
This function should be called after
solve()has finished.The value
nconvmay be different from the number of requested solutionsnev, but not larger thanncv, seesetDimensions().See also
- getConvergedReason()#
Get the reason why the
solve()iteration was stopped.Not collective.
- Returns:
Negative value indicates diverged, positive value converged.
- Return type:
See also
- getConvergenceTest()#
Get how to compute the error estimate used in the convergence test.
Not collective.
- Returns:
The method used to compute the error estimate used in the convergence test.
- Return type:
See also
- getDS()#
Get the direct solver associated to the eigensolver.
Not collective.
- Returns:
The direct solver context.
- Return type:
- getDimensions()#
Get number of eigenvalues to compute and the dimension of the subspace.
Not collective.
- Returns:
- Return type:
See also
- getEigenpair(i, Vr=None, Vi=None)#
Get the i-th solution of the eigenproblem as computed by
solve().Collective.
The solution consists of both the eigenvalue and the eigenvector.
- Parameters:
- Returns:
e – The computed eigenvalue. It will be a real variable in case of a Hermitian or generalized Hermitian eigenproblem. Otherwise it will be a complex variable (possibly with zero imaginary part).
- Return type:
Notes
The index
ishould be a value between0andnconv-1(seegetConverged()). Eigenpairs are indexed according to the ordering criterion established withsetWhichEigenpairs().The 2-norm of the eigenvector is one unless the problem is generalized Hermitian. In this case the eigenvector is normalized with respect to the norm defined by the B matrix.
See also
- getEigenvalue(i)#
Get the i-th eigenvalue as computed by
solve().Not collective.
- Parameters:
i (int) – Index of the solution to be obtained.
- Returns:
The computed eigenvalue. It will be a real variable in case of a Hermitian or generalized Hermitian eigenproblem. Otherwise it will be a complex variable (possibly with zero imaginary part).
- Return type:
Notes
The index
ishould be a value between0andnconv-1(seegetConverged()). Eigenpairs are indexed according to the ordering criterion established withsetWhichEigenpairs().See also
getConverged,setWhichEigenpairs,getEigenpair,EPSGetEigenvalue
- getEigenvalueComparison()#
Get the eigenvalue comparison function.
Not collective.
- Returns:
The eigenvalue comparison function.
- Return type:
See also
- getEigenvector(i, Vr=None, Vi=None)#
Get the i-th right eigenvector as computed by
solve().Collective.
- Parameters:
- Return type:
Notes
The index
ishould be a value between0andnconv-1(seegetConverged()). Eigenpairs are indexed according to the ordering criterion established withsetWhichEigenpairs().The 2-norm of the eigenvector is one unless the problem is generalized Hermitian. In this case the eigenvector is normalized with respect to the norm defined by the B matrix.
See also
getConverged,setWhichEigenpairs,getEigenpair,EPSGetEigenvector
- getErrorEstimate(i)#
Get the error estimate associated to the i-th computed eigenpair.
Not collective.
- Parameters:
i (int) – Index of the solution to be considered.
- Returns:
Error estimate.
- Return type:
Notes
This is the error estimate used internally by the eigensolver. The actual error bound can be computed with
computeError().See also
- getExtraction()#
Get the extraction type used by the EPS object.
Not collective.
- Returns:
The method of extraction.
- Return type:
See also
- getGDBOrth()#
Get the orthogonalization used in the search subspace.
Not collective.
Get the orthogonalization used in the search subspace in case of generalized Hermitian problems.
- Returns:
Whether to B-orthogonalize the search subspace.
- Return type:
See also
- getGDBlockSize()#
Get the number of vectors to be added to the searching space.
Not collective.
Get the number of vectors to be added to the searching space in every iteration.
- Returns:
The number of vectors added to the search space in every iteration.
- Return type:
See also
- getGDDoubleExpansion()#
Get a flag indicating whether the double expansion variant is active.
Not collective.
Get a flag indicating whether the double expansion variant has been activated or not.
- Returns:
Trueif using double expansion.- Return type:
See also
- getGDInitialSize()#
Get the initial size of the searching space.
Not collective.
- Returns:
The number of vectors of the initial searching subspace.
- Return type:
See also
- getGDKrylovStart()#
Get a flag indicating if the search subspace is started with a Krylov basis.
Not collective.
- Returns:
Trueif starting the search subspace with a Krylov basis.- Return type:
See also
- getGDRestart()#
Get the number of vectors of the search space after restart.
Not collective.
Get the number of vectors of the search space after restart and the number of vectors saved from the previous iteration.
- Returns:
- Return type:
See also
- getInterval()#
Get the computational interval for spectrum slicing.
Not collective.
- Returns:
- Return type:
Notes
If the interval was not set by the user, then zeros are returned.
See also
- getInvariantSubspace()#
Get an orthonormal basis of the computed invariant subspace.
Collective.
- Returns:
Basis of the invariant subspace.
- Return type:
Notes
This function should be called after
solve()has finished.The returned vectors span an invariant subspace associated with the computed eigenvalues. An invariant subspace \(X\) of \(A\) satisfies \(A x \in X\), for all \(x \in X\) (a similar definition applies for generalized eigenproblems).
See also
- getIterationNumber()#
Get the current iteration number.
Not collective.
If the call to
solve()is complete, then it returns the number of iterations carried out by the solution method.- Returns:
Iteration number.
- Return type:
See also
- getJDBOrth()#
Get the orthogonalization used in the search subspace.
Not collective.
Get the orthogonalization used in the search subspace in case of generalized Hermitian problems.
- Returns:
Whether to B-orthogonalize the search subspace.
- Return type:
See also
- getJDBlockSize()#
Get the number of vectors to be added to the searching space.
Not collective.
Get the number of vectors to be added to the searching space in every iteration.
- Returns:
The number of vectors added to the search space in every iteration.
- Return type:
See also
- getJDConstCorrectionTol()#
Get the flag indicating if the dynamic stopping is being used.
Not collective.
- Returns:
Trueif the dynamic stopping criterion is not being used.- Return type:
- getJDFix()#
Get the threshold for changing the target in the correction equation.
Not collective.
- Returns:
The threshold for changing the target.
- Return type:
See also
- getJDInitialSize()#
Get the initial size of the searching space.
Not collective.
- Returns:
The number of vectors of the initial searching subspace.
- Return type:
See also
- getJDKrylovStart()#
Get a flag indicating if the search subspace is started with a Krylov basis.
Not collective.
- Returns:
Trueif starting the search subspace with a Krylov basis.- Return type:
See also
- getJDRestart()#
Get the number of vectors of the search space after restart.
Not collective.
Get the number of vectors of the search space after restart and the number of vectors saved from the previous iteration.
- Returns:
- Return type:
See also
- getKrylovSchurBSEType()#
Get the method used for BSE structured eigenproblems (Krylov-Schur).
Not collective.
- Returns:
The BSE method.
- Return type:
- getKrylovSchurDetectZeros()#
Get the flag that enforces zero detection in spectrum slicing.
Not collective.
- Returns:
The zero detection flag.
- Return type:
- getKrylovSchurDimensions()#
Get the dimensions used for each subsolve step (spectrum slicing).
Not collective.
- Returns:
- Return type:
- getKrylovSchurInertias()#
Get the values of the shifts and their corresponding inertias.
Not collective.
Get the values of the shifts and their corresponding inertias in case of doing spectrum slicing for a computational interval.
- Returns:
- Return type:
Notes
This call makes sense only for spectrum slicing runs, that is, when an interval has been given with
setInterval()andSINVERTis set.If called after
solve(), all shifts used internally by the solver are returned (including both endpoints and any intermediate ones). If called beforesolve()and aftersetUp()then only the information of the endpoints of subintervals is available.
- getKrylovSchurKSP()#
Get the linear solver object associated with the internal
EPSobject.Collective.
Get the linear solver object associated with the internal
EPSobject in case of doing spectrum slicing for a computational interval.- Returns:
The linear solver object.
- Return type:
Notes
This call makes sense only for spectrum slicing runs, that is, when an interval has been given with
setInterval()andSINVERTis set.When invoked to compute all eigenvalues in an interval with spectrum slicing,
KRYLOVSCHURcreates anotherEPSobject internally that is used to compute eigenvalues by chunks near selected shifts. This function allows access to theKSPobject associated to this internalEPSobject.In case of having more than one partition, the returned
KSPwill be different in MPI processes belonging to different partitions. Hence, if required,setKrylovSchurPartitions()must be called BEFORE this function.
- getKrylovSchurLocking()#
Get the locking flag used in the Krylov-Schur method.
Not collective.
- Returns:
The locking flag.
- Return type:
- getKrylovSchurPartitions()#
Get the number of partitions of the communicator (spectrum slicing).
Not collective.
- Returns:
The number of partitions.
- Return type:
- getKrylovSchurRestart()#
Get the restart parameter used in the Krylov-Schur method.
Not collective.
- Returns:
The number of vectors to be kept at restart.
- Return type:
- getKrylovSchurSubcommInfo()#
Get information related to the case of doing spectrum slicing.
Collective on the subcommunicator.
Get information related to the case of doing spectrum slicing for a computational interval with multiple communicators.
- Returns:
k (
int) – Index of the subinterval for the calling process.n (
int) – Number of eigenvalues found in thek-th subinterval.v (
petsc4py.PETSc.Vec) – A vector owned by processes in the subcommunicator with dimensions compatible for locally computed eigenvectors.
- Return type:
Notes
This call makes sense only for spectrum slicing runs, that is, when an interval has been given with
setInterval()andSINVERTis set.
- getKrylovSchurSubcommMats()#
Get the eigenproblem matrices stored in the subcommunicator.
Collective on the subcommunicator.
Get the eigenproblem matrices stored internally in the subcommunicator to which the calling process belongs.
- Returns:
A (
petsc4py.PETSc.Mat) – The matrix associated with the eigensystem.B (
petsc4py.PETSc.Mat) – The second matrix in the case of generalized eigenproblems.
- Return type:
Notes
This call makes sense only for spectrum slicing runs, that is, when an interval has been given with
setInterval()andSINVERTis set. And is relevant only when the number of partitions (setKrylovSchurPartitions()) is larger than one.This is the analog of
getOperators(), but returns the matrices distributed differently (in the subcommunicator rather than in the parent communicator).These matrices should not be modified by the user.
- getKrylovSchurSubcommPairs(i, v=None)#
Get the i-th eigenpair stored in the multi-communicator of the process.
Collective on the subcommunicator (if v is given).
Get the i-th eigenpair stored internally in the multi-communicator to which the calling process belongs.
- Parameters:
- Returns:
The computed eigenvalue.
- Return type:
Notes
This call makes sense only for spectrum slicing runs, that is, when an interval has been given with
setInterval()andSINVERTis set. And is relevant only when the number of partitions (setKrylovSchurPartitions()) is larger than one.Argument
vmust be a validVecobject, created by callinggetKrylovSchurSubcommInfo().The index
ishould be a value between0andn-1, wherenis the number of vectors in the local subinterval, seegetKrylovSchurSubcommInfo().
- getKrylovSchurSubintervals()#
Get the points that delimit the subintervals.
Not collective.
Get the points that delimit the subintervals used in spectrum slicing with several partitions.
- Returns:
Real values specifying subintervals.
- Return type:
Notes
This call makes sense only for spectrum slicing runs, that is, when an interval has been given with
setInterval()andSINVERTis set.If the user passed values with
setKrylovSchurSubintervals(), then the same values are returned here. Otherwise, the values computed internally are obtained.
- getLOBPCGBlockSize()#
Get the block size used in the LOBPCG method.
Not collective.
- Returns:
The block size.
- Return type:
See also
- getLOBPCGLocking()#
Get the locking flag used in the LOBPCG method.
Not collective.
- Returns:
The locking flag.
- Return type:
See also
- getLOBPCGRestart()#
Get the restart parameter used in the LOBPCG method.
Not collective.
- Returns:
The restart parameter.
- Return type:
See also
- getLanczosReorthogType()#
Get the type of reorthogonalization used during the Lanczos iteration.
Not collective.
- Returns:
The type of reorthogonalization.
- Return type:
See also
- getLeftEigenvector(i, Wr=None, Wi=None)#
Get the i-th left eigenvector as computed by
solve().Collective.
- Parameters:
- Return type:
Notes
The index
ishould be a value between0andnconv-1(seegetConverged()). Eigensolutions are indexed according to the ordering criterion established withsetWhichEigenpairs().Left eigenvectors are available only if the twosided flag was set with
setTwoSided().
- getLyapIIRanks()#
Get the rank values used for the Lyapunov step.
Not collective.
- Returns:
- Return type:
See also
- getMonitor()#
Get the list of monitor functions.
Not collective.
- Returns:
The list of monitor functions.
- Return type:
See also
- getOperators()#
Get the matrices associated with the eigenvalue problem.
Collective.
- Returns:
A (
petsc4py.PETSc.Mat) – The matrix associated with the eigensystem.B (
petsc4py.PETSc.Mat) – The second matrix in the case of generalized eigenproblems.
- Return type:
See also
- getOptionsPrefix()#
Get the prefix used for searching for all EPS options in the database.
Not collective.
- Returns:
The prefix string set for this EPS object.
- Return type:
- getPowerShiftType()#
Get the type of shifts used during the power iteration.
Not collective.
- Returns:
The type of shift.
- Return type:
See also
- getProblemType()#
Get the problem type from the EPS object.
Not collective.
- Returns:
The problem type that was previously set.
- Return type:
See also
- getPurify()#
Get the flag indicating whether purification is activated or not.
Not collective.
- Returns:
Whether purification is activated or not.
- Return type:
See also
- getRG()#
Get the region object associated to the eigensolver.
Not collective.
- Returns:
The region context.
- Return type:
- getRQCGReset()#
Get the reset parameter used in the RQCG method.
Not collective.
- Returns:
The number of iterations between resets.
- Return type:
See also
- getST()#
Get the spectral transformation object associated to the eigensolver.
Not collective.
- Returns:
The spectral transformation.
- Return type:
- getStoppingTest()#
Get the stopping test function.
Not collective.
- Returns:
The stopping test function.
- Return type:
See also
- getTarget()#
Get the value of the target.
Not collective.
- Returns:
The value of the target.
- Return type:
Notes
If the target was not set by the user, then zero is returned.
See also
- getThreshold()#
Get the threshold used in the threshold stopping test.
Not collective.
- Returns:
- Return type:
See also
- getTolerances()#
Get the tolerance and max. iter. count used for convergence tests.
Not collective.
Get the tolerance and iteration limit used by the default EPS convergence tests.
- Returns:
- Return type:
See also
- getTrackAll()#
Get the flag indicating if all residual norms must be computed or not.
Not collective.
- Returns:
Whether the solver computes all residuals or not.
- Return type:
See also
- getTrueResidual()#
Get the flag indicating if true residual must be computed explicitly.
Not collective.
- Returns:
Whether the solver computes true residuals or not.
- Return type:
See also
- getTwoSided()#
Get the flag indicating if a two-sided variant of the algorithm is being used.
Not collective.
- Returns:
Whether the two-sided variant is to be used or not.
- Return type:
See also
- getType()#
Get the EPS type of this object.
Not collective.
- Returns:
The solver currently being used.
- Return type:
See also
- getWhichEigenpairs()#
Get which portion of the spectrum is to be sought.
Not collective.
- Returns:
The portion of the spectrum to be sought by the solver.
- Return type:
See also
- isGeneralized()#
Tell if the EPS object corresponds to a generalized eigenproblem.
Not collective.
- Returns:
Trueif the problem is generalized.- Return type:
See also
- isHermitian()#
Tell if the EPS object corresponds to a Hermitian eigenproblem.
Not collective.
- Returns:
Trueif the problem is Hermitian.- Return type:
See also
- isPositive()#
Eigenproblem requiring a positive (semi-) definite matrix \(B\).
Not collective.
Tell if the EPS corresponds to an eigenproblem requiring a positive (semi-) definite matrix \(B\).
- Returns:
Trueif the problem is positive (semi-) definite.- Return type:
See also
- isStructured()#
Tell if the EPS object corresponds to a structured eigenvalue problem.
Not collective.
- Returns:
Trueif the problem is structured.- Return type:
Notes
The result will be
Trueif the problem type has been set to some structured type such asBSE. This is independent of whether the input matrix has been built with a certain structure with a helper function.See also
- reset()#
Reset the EPS object.
Collective.
See also
Source code at slepc4py/SLEPc/EPS.pyx:384
- Return type:
- setArbitrarySelection(arbitrary, args=None, kargs=None)#
Set an arbitrary selection criterion function.
Logically collective.
Set a function to look for eigenvalues according to an arbitrary selection criterion. This criterion can be based on a computation involving the current eigenvector approximation.
- setArnoldiDelayed(delayed)#
Set (toggle) delayed reorthogonalization in the Arnoldi iteration.
Logically collective.
Notes
This call is only relevant if the type was set to
EPS.Type.ARNOLDIwithsetType().Delayed reorthogonalization is an aggressive optimization for the Arnoldi eigensolver than may provide better scalability, but sometimes makes the solver converge more slowly compared to the default algorithm.
See also
- setBV(bv)#
Set a basis vectors object associated to the eigensolver.
Collective.
- setBalance(balance=None, iterations=None, cutoff=None)#
Set the balancing technique to be used by the eigensolver.
Logically collective.
- Parameters:
- Return type:
Notes
When balancing is enabled, the solver works implicitly with matrix \(DAD^{-1}\), where \(D\) is an appropriate diagonal matrix. This improves the accuracy of the computed results in some cases.
Balancing makes sense only for non-Hermitian problems when the required precision is high (i.e., with a small tolerance).
By default, balancing is disabled. The two-sided method is much more effective than the one-sided counterpart, but it requires the system matrices to have the
Mat.multTranspose()operation defined.The parameter
iterationsis the number of iterations performed by the method. Thecutoffvalue is used only in the two-side variant.See also
- setCISSExtraction(extraction)#
Set the extraction technique used in the CISS solver.
Logically collective.
- Parameters:
extraction (CISSExtraction) – The extraction technique.
- Return type:
See also
- setCISSQuadRule(quad)#
Set the quadrature rule used in the CISS solver.
Logically collective.
- Parameters:
quad (CISSQuadRule) – The quadrature rule.
- Return type:
See also
- setCISSRefinement(inner=None, blsize=None)#
Set the values of various refinement parameters in the CISS solver.
Logically collective.
- Parameters:
- Return type:
See also
- setCISSSizes(ip=None, bs=None, ms=None, npart=None, bsmax=None, realmats=False)#
Set the values of various size parameters in the CISS solver.
Logically collective.
- Parameters:
- Return type:
Notes
The default number of partitions is 1. This means the internal
petsc4py.PETSc.KSPobject is shared among all processes of theEPScommunicator. Otherwise, the communicator is split intonpartcommunicators, so thatnpartpetsc4py.PETSc.KSPsolves proceed simultaneously.See also
getCISSSizes,setCISSThreshold,setCISSRefinement,EPSCISSSetSizes
- setCISSThreshold(delta=None, spur=None)#
Set the values of various threshold parameters in the CISS solver.
Logically collective.
- Parameters:
- Return type:
See also
- setCISSUseST(usest)#
Set a flag indicating that the CISS solver will use the
STobject.Logically collective.
Notes
When this option is set, the linear solves can be configured by setting options for the
petsc4py.PETSc.KSPobject obtained withST.getKSP(). Otherwise, severalpetsc4py.PETSc.KSPobjects are created, which can be accessed withgetCISSKSPs().The default is to use the
ST, unless several partitions have been specified, seesetCISSSizes().See also
- setConvergenceTest(conv)#
Set how to compute the error estimate used in the convergence test.
Logically collective.
- Parameters:
conv (Conv) – The method used to compute the error estimate used in the convergence test.
- Return type:
See also
- setDS(ds)#
Set a direct solver object associated to the eigensolver.
Collective.
- setDeflationSpace(space)#
Set vectors to form a basis of the deflation space.
Collective.
Notes
When a deflation space is given, the eigensolver seeks the eigensolution in the restriction of the problem to the orthogonal complement of this space. This can be used for instance in the case that an invariant subspace is known beforehand (such as the nullspace of the matrix).
These vectors do not persist from one
solve()call to the other, so the deflation space should be set every time.The vectors do not need to be mutually orthonormal, since they are explicitly orthonormalized internally.
See also
- setDimensions(nev=None, ncv=None, mpd=None)#
Set number of eigenvalues to compute and the dimension of the subspace.
Logically collective.
- Parameters:
- Return type:
Notes
Use
DETERMINEforncvandmpdto assign a reasonably good value, which is dependent on the solution method.The parameters
ncvandmpdare intimately related, so that the user is advised to set one of them at most. Normal usage is the following:In cases where
nevis small, the user setsncv(a reasonable default is 2 *nev).In cases where
nevis large, the user setsmpd.
The value of
ncvshould always be betweennevand (nev+mpd), typicallyncv=nev+mpd. Ifnevis not too large,mpd=nevis a reasonable choice, otherwise a smaller value should be used.When computing all eigenvalues in an interval, see
setInterval(), these parameters lose relevance, and tuning must be done withsetKrylovSchurDimensions().
- setEigenvalueComparison(comparison, args=None, kargs=None)#
Set an eigenvalue comparison function.
Logically collective.
Notes
This eigenvalue comparison function is used when
setWhichEigenpairs()is set toEPS.Which.USER.
- setExtraction(extraction)#
Set the extraction type used by the eigensolver.
Logically collective.
- Parameters:
extraction (Extraction) – The extraction method to be used by the solver.
- Return type:
Notes
Not all eigensolvers support all types of extraction.
By default, a standard Rayleigh-Ritz extraction is used. Other extractions may be useful when computing interior eigenvalues.
Harmonic-type extractions are used in combination with a target, see
setTarget().See also
- setFromOptions()#
Set EPS options from the options database.
Collective.
Notes
To see all options, run your program with the
-helpoption.This routine must be called before
setUp()if the user is to be allowed to set the solver type.See also
Source code at slepc4py/SLEPc/EPS.pyx:532
- Return type:
- setGDBOrth(borth)#
Set the orthogonalization that will be used in the search subspace.
Logically collective.
Set the orthogonalization that will be used in the search subspace in case of generalized Hermitian problems.
See also
- setGDBlockSize(bs)#
Set the number of vectors to be added to the searching space.
Logically collective.
Set the number of vectors to be added to the searching space in every iteration.
- Parameters:
bs (int) – The number of vectors added to the search space in every iteration.
- Return type:
See also
- setGDDoubleExpansion(doubleexp)#
Set that the search subspace is expanded with double expansion.
Logically collective.
Notes
In the double expansion variant the search subspace is expanded with \(K [A x, B x]\) (double expansion) instead of the classic \(K r\), where \(K\) is the preconditioner, \(x\) the selected approximate eigenvector and \(r\) its associated residual vector.
See also
- setGDInitialSize(initialsize)#
Set the initial size of the searching space.
Logically collective.
- Parameters:
initialsize (int) – The number of vectors of the initial searching subspace.
- Return type:
Notes
If the flag in
setGDKrylovStart()is set toFalseand the user provides vectors withsetInitialSpace(), up toinitialsizevectors will be used; and if the provided vectors are not enough, the solver completes the subspace with random vectors. In case thesetGDKrylovStart()flag isTrue, the solver gets the first vector provided by the user or, if not available, a random vector, and expands the Krylov basis up toinitialsizevectors.See also
- setGDKrylovStart(krylovstart=True)#
Set (toggle) starting the search subspace with a Krylov basis.
Logically collective.
- Parameters:
krylovstart (bool) –
Trueif starting the search subspace with a Krylov basis.- Return type:
See also
- setGDRestart(minv=None, plusk=None)#
Set the number of vectors of the search space after restart.
Logically collective.
Set the number of vectors of the search space after restart and the number of vectors saved from the previous iteration.
- Parameters:
- Return type:
See also
- setInitialSpace(space)#
Set the initial space from which the eigensolver starts to iterate.
Collective.
Notes
Some solvers start to iterate on a single vector (initial vector). In that case, only the first vector is taken into account and the other vectors are ignored. But other solvers such as
SUBSPACEare able to make use of the whole initial subspace as an initial guess.These vectors do not persist from one
solve()call to the other, so the initial space should be set every time.The vectors do not need to be mutually orthonormal, since they are explicitly orthonormalized internally.
Common usage of this function is when the user can provide a rough approximation of the wanted eigenspace. Then, convergence may be faster.
- setInterval(inta, intb)#
Set the computational interval for spectrum slicing.
Logically collective.
- Parameters:
- Return type:
Notes
Spectrum slicing is a technique employed for computing all eigenvalues of symmetric eigenproblems in a given interval. This function provides the interval to be considered. It must be used in combination with
EPS.Which.ALL, seesetWhichEigenpairs().A computational interval is also needed when using polynomial filters, see
STFILTER.See also
- setJDBOrth(borth)#
Set the orthogonalization that will be used in the search subspace.
Logically collective.
Set the orthogonalization that will be used in the search subspace in case of generalized Hermitian problems.
See also
- setJDBlockSize(bs)#
Set the number of vectors to be added to the searching space.
Logically collective.
Set the number of vectors to be added to the searching space in every iteration.
- Parameters:
bs (int) – The number of vectors added to the search space in every iteration.
- Return type:
See also
- setJDConstCorrectionTol(constant)#
Deactivate the dynamic stopping criterion.
Logically collective.
- Parameters:
constant (bool) – If
False, thepetsc4py.PETSc.KSPrelative tolerance is set to0.5**i.- Return type:
Notes
If this flag is set to
False, then thepetsc4py.PETSc.KSPrelative tolerance is dynamically set to0.5**i, whereiis the number ofEPSiterations since the last converged value. By the default, a constant tolerance is used.
- setJDFix(fix)#
Set the threshold for changing the target in the correction equation.
Logically collective.
Notes
The target in the correction equation is fixed at the first iterations. When the norm of the residual vector is lower than the
fixvalue, the target is set to the corresponding eigenvalue.See also
- setJDInitialSize(initialsize)#
Set the initial size of the searching space.
Logically collective.
- Parameters:
initialsize (int) – The number of vectors of the initial searching subspace.
- Return type:
Notes
If the flag in
setJDKrylovStart()is set toFalseand the user provides vectors withsetInitialSpace(), up toinitialsizevectors will be used; and if the provided vectors are not enough, the solver completes the subspace with random vectors. In case thesetJDKrylovStart()flag isTrue, the solver gets the first vector provided by the user or, if not available, a random vector, and expands the Krylov basis up toinitialsizevectors.See also
- setJDKrylovStart(krylovstart=True)#
Set (toggle) starting the search subspace with a Krylov basis.
Logically collective.
- Parameters:
krylovstart (bool) –
Trueif starting the search subspace with a Krylov basis.- Return type:
See also
- setJDRestart(minv=None, plusk=None)#
Set the number of vectors of the search space after restart.
Logically collective.
Set the number of vectors of the search space after restart and the number of vectors saved from the previous iteration.
- Parameters:
- Return type:
See also
- setKrylovSchurBSEType(bse)#
Set the Krylov-Schur variant used for BSE structured eigenproblems.
Logically collective.
- Parameters:
bse (KrylovSchurBSEType) – The BSE method.
- Return type:
Notes
This call is only relevant if the type was set to
EPS.Type.KRYLOVSCHURwithsetType()and the problem type toEPS.ProblemType.BSEwithsetProblemType().
- setKrylovSchurDetectZeros(detect)#
Set the flag that enforces zero detection in spectrum slicing.
Logically collective.
Set a flag to enforce the detection of zeros during the factorizations throughout the spectrum slicing computation.
Notes
This call makes sense only for spectrum slicing runs, that is, when an interval has been given with
setInterval()andSINVERTis set.A zero in the factorization indicates that a shift coincides with an eigenvalue.
This flag is turned off by default, and may be necessary in some cases, especially when several partitions are being used. This feature currently requires an external package for factorizations with support for zero detection, e.g., MUMPS.
- setKrylovSchurDimensions(nev=None, ncv=None, mpd=None)#
Set the dimensions used for each subsolve step (spectrum slicing).
Logically collective.
- Parameters:
- Return type:
Notes
This call makes sense only for spectrum slicing runs, that is, when an interval has been given with
setInterval()andSINVERTis set.The meaning of the parameters is the same as in
setDimensions(), but the ones here apply to every subsolve done by the childEPSobject.
- setKrylovSchurLocking(lock)#
Set (toggle) locking/non-locking variants of the Krylov-Schur method.
Logically collective.
Notes
The default is to lock converged eigenpairs when the method restarts. This behavior can be changed so that all directions are kept in the working subspace even if already converged to working accuracy (the non-locking variant).
- setKrylovSchurPartitions(npart)#
Set the number of partitions of the communicator (spectrum slicing).
Logically collective.
Set the number of partitions for the case of doing spectrum slicing for a computational interval with the communicator split in several sub-communicators.
Notes
This call makes sense only for spectrum slicing runs, that is, when an interval has been given with
setInterval()andSINVERTis set.By default,
npart=1so all processes in the communicator participate in the processing of the whole interval. Ifnpart>1then the interval is divided intonpartsubintervals, each of them being processed by a subset of processes.The interval is split proportionally unless the separation points are specified with
setKrylovSchurSubintervals().
- setKrylovSchurRestart(keep)#
Set the restart parameter for the Krylov-Schur method.
Logically collective.
It is the proportion of basis vectors that must be kept after restart.
Notes
Allowed values are in the range [0.1,0.9]. The default is 0.5.
- setKrylovSchurSubintervals(subint)#
Set the subinterval boundaries.
Logically collective.
Set the subinterval boundaries for spectrum slicing with a computational interval with several partitions.
Notes
This call makes sense only for spectrum slicing runs, that is, when an interval has been given with
setInterval()andSINVERTis set.This function must be called after
setKrylovSchurPartitions(). Fornpartpartitions, the argumentsubintmust containnpart+1real values sorted in ascending order:subint_0,subint_1, …,subint_npart, where the first and last values must coincide with the interval endpoints set withsetInterval(). The subintervals are then defined by two consecutive points:[subint_0,subint_1],[subint_1,subint_2], and so on.
- setLOBPCGBlockSize(bs)#
Set the block size of the LOBPCG method.
Logically collective.
See also
- setLOBPCGLocking(lock)#
Toggle between locking and non-locking (LOBPCG method).
Logically collective.
Notes
This flag refers to soft locking (converged vectors within the current block iterate), since hard locking is always used (when
nevis larger than the block size).See also
- setLOBPCGRestart(restart)#
Set the restart parameter for the LOBPCG method.
Logically collective.
- Parameters:
restart (float) – The percentage of the block of vectors to force a restart.
- Return type:
Notes
The meaning of this parameter is the proportion of vectors within the current block iterate that must have converged in order to force a restart with hard locking. Allowed values are in the range [0.1,1.0]. The default is 0.9.
See also
- setLanczosReorthogType(reorthog)#
Set the type of reorthogonalization used during the Lanczos iteration.
Logically collective.
- Parameters:
reorthog (LanczosReorthogType) – The type of reorthogonalization.
- Return type:
Notes
This call is only relevant if the type was set to
EPS.Type.LANCZOSwithsetType().See also
- setLeftInitialSpace(space)#
Set a left initial space from which the eigensolver starts to iterate.
Collective.
Notes
Left initial vectors are used to initiate the left search space in two-sided eigensolvers. Users should pass here an approximation of the left eigenspace, if available.
The same comments in
setInitialSpace()are applicable here.See also
- setLyapIIRanks(rkc=None, rkl=None)#
Set the ranks used in the solution of the Lyapunov equation.
Logically collective.
- Parameters:
- Return type:
Notes
Lyapunov inverse iteration needs to solve a large-scale Lyapunov equation at each iteration of the eigensolver. For this, an iterative solver (
LME) is used, which requires to prescribe the rank of the solution matrix \(X\). This is the meaning of parameterrkl. Later, this matrix is compressed into another matrix of rankrkc. If not provided,rklis a small multiple ofrkc.See also
- setMonitor(monitor, args=None, kargs=None)#
Append a monitor function to the list of monitors.
Logically collective.
See also
- setOperators(A, B=None)#
Set the matrices associated with the eigenvalue problem.
Collective.
- Parameters:
- Return type:
Notes
It must be called before
setUp(). If it is called again aftersetUp()and the matrix sizes have changed then theEPSobject is reset.For structured eigenproblem types such as
BSE, seesetProblemType(), the provided matrices must have been created with the corresponding helper function, i.e.,createMatBSE().See also
getOperators,solve,setUp,reset,setProblemType,EPSSetOperators
- setOptionsPrefix(prefix=None)#
Set the prefix used for searching for all EPS options in the database.
Logically collective.
- Parameters:
prefix (str | None) – The prefix string to prepend to all EPS option requests.
- Return type:
Notes
A hyphen (-) must NOT be given at the beginning of the prefix name. The first character of all runtime options is AUTOMATICALLY the hyphen.
For example, to distinguish between the runtime options for two different EPS contexts, one could call:
E1.setOptionsPrefix("eig1_") E2.setOptionsPrefix("eig2_")
- setPowerShiftType(shift)#
Set the type of shifts used during the power iteration.
Logically collective.
This can be used to emulate the Rayleigh Quotient Iteration (RQI) method.
- Parameters:
shift (PowerShiftType) – The type of shift.
- Return type:
Notes
This call is only relevant if the type was set to
EPS.Type.POWERwithsetType().By default, shifts are constant (
EPS.PowerShiftType.CONSTANT) and the iteration is the simple power method (or inverse iteration if a shift-and-invert transformation is being used).A variable shift can be specified (
EPS.PowerShiftType.RAYLEIGHorEPS.PowerShiftType.WILKINSON). In this case, the iteration behaves rather like a cubic converging method as RQI.See also
- setProblemType(problem_type)#
Set the type of the eigenvalue problem.
Logically collective.
- Parameters:
problem_type (ProblemType) – The problem type to be set.
- Return type:
Notes
This function must be used to instruct SLEPc to exploit symmetry or other kind of structure. If no problem type is specified, by default a non-Hermitian problem is assumed (either standard or generalized). If the user knows that the problem is Hermitian (i.e., \(A=A^*\)) or generalized Hermitian (i.e., \(A=A^*\), \(B=B^*\), and \(B\) positive definite) then it is recommended to set the problem type so that eigensolver can exploit these properties.
If the user does not call this function, the solver will use a reasonable guess.
For structured problem types such as
BSE, the matrices passed in viasetOperators()must have been created with the corresponding helper function, i.e.,createMatBSE().See also
setOperators,createMatBSE,getProblemType,EPSSetProblemType
- setPurify(purify=True)#
Set (toggle) eigenvector purification.
Logically collective.
Notes
By default, eigenvectors of generalized symmetric eigenproblems are purified in order to purge directions in the nullspace of matrix \(B\). If the user knows that \(B\) is non-singular, then purification can be safely deactivated and some computational cost is avoided (this is particularly important in interval computations).
See also
- setRG(rg)#
Set a region object associated to the eigensolver.
Collective.
- setRQCGReset(nrest)#
Set the reset parameter of the RQCG iteration.
Logically collective.
Notes
Every
nrestiterations the solver performs a Rayleigh-Ritz projection step.See also
- setST(st)#
Set a spectral transformation object associated to the eigensolver.
Collective.
- setStoppingTest(stopping, args=None, kargs=None)#
Set a function to decide when to stop the outer iteration of the eigensolver.
Logically collective.
See also
- setTarget(target)#
Set the value of the target.
Logically collective.
Notes
The target is a scalar value used to determine the portion of the spectrum of interest. It is used in combination with
setWhichEigenpairs().When PETSc is built with real scalars, it is not possible to specify a complex target.
See also
- setThreshold(thres, rel=False)#
Set the threshold used in the threshold stopping test.
Logically collective.
- Parameters:
- Return type:
Notes
This function internally sets a special stopping test based on the threshold, where eigenvalues are computed in sequence until one of the computed eigenvalues is below/above the threshold (depending on whether largest or smallest eigenvalues are computed). The details are given in
EPSSetThreshold.See also
- setTolerances(tol=None, max_it=None)#
Set the tolerance and max. iter. used by the default EPS convergence tests.
Logically collective.
- Parameters:
- Return type:
Notes
Use
DETERMINEformax_itto assign a reasonably good value, which is dependent on the solution method.See also
- setTrackAll(trackall)#
Set if the solver must compute the residual of all approximate eigenpairs.
Logically collective.
See also
- setTrueResidual(trueres)#
Set if the solver must compute the true residual explicitly or not.
Logically collective.
See also
- setTwoSided(twosided)#
Set to use a two-sided variant that also computes left eigenvectors.
Logically collective.
Notes
If the user sets
twosidedtoTruethen the solver uses a variant of the algorithm that computes both right and left eigenvectors. This is usually much more costly. This option is not available in all solvers.When using two-sided solvers, the problem matrices must have both the
Mat.multandMat.multTransposeoperations defined.See also
- setType(eps_type)#
Set the particular solver to be used in the EPS object.
Logically collective.
Notes
The default is
KRYLOVSCHUR. Normally, it is best to usesetFromOptions()and then set the EPS type from the options database rather than by using this routine. Using the options database provides the user with maximum flexibility in evaluating the different available methods.See also
- setUp()#
Set up all the internal data structures.
Collective.
Notes
Sets up all the internal data structures necessary for the execution of the eigensolver. This includes the setup of the internal
STobject.This function need not be called explicitly in most cases, since
solve()calls it. It can be useful when one wants to measure the set-up time separately from the solve time.See also
Source code at slepc4py/SLEPc/EPS.pyx:1875
- Return type:
- setWhichEigenpairs(which)#
Set which portion of the spectrum is to be sought.
Logically collective.
- Parameters:
which (Which) – The portion of the spectrum to be sought by the solver.
- Return type:
Notes
Not all eigensolvers implemented in EPS account for all the possible values. Also, some values make sense only for certain types of problems. If SLEPc is compiled for real numbers
EPS.Which.LARGEST_IMAGINARYandEPS.Which.SMALLEST_IMAGINARYuse the absolute value of the imaginary part for eigenvalue selection.The target is a scalar value provided with
setTarget().The criterion
EPS.Which.TARGET_IMAGINARYis available only in case PETSc and SLEPc have been built with complex scalars.EPS.Which.ALLis intended for use in combination with an interval (seesetInterval()), when all eigenvalues within the interval are requested, or in the context of theEPS.Type.CISSsolver for computing all eigenvalues in a region.See also
setTarget,setInterval,getWhichEigenpairs,EPSSetWhichEigenpairs
- solve()#
Solve the eigensystem.
Collective.
Notes
The problem matrices are specified with
setOperators().solve()will return without generating an error regardless of whether all requested solutions were computed or not. CallgetConverged()to get the actual number of computed solutions, andgetConvergedReason()to determine if the solver converged or failed and why.See also
setUp,setOperators,getConverged,getConvergedReason,EPSSolveSource code at slepc4py/SLEPc/EPS.pyx:1897
- Return type:
- updateKrylovSchurSubcommMats(s=1.0, a=1.0, Au=None, t=1.0, b=1.0, Bu=None, structure=None, globalup=False)#
Update the eigenproblem matrices stored internally in the communicator.
Collective.
Update the eigenproblem matrices stored internally in the subcommunicator to which the calling process belongs.
- Parameters:
s (Scalar) – Scalar that multiplies the existing A matrix.
a (Scalar) – Scalar used in the axpy operation on A.
Au (petsc4py.PETSc.Mat | None) – The matrix used in the axpy operation on A.
t (Scalar) – Scalar that multiplies the existing B matrix.
b (Scalar) – Scalar used in the axpy operation on B.
Bu (petsc4py.PETSc.Mat | None) – The matrix used in the axpy operation on B.
structure (petsc4py.PETSc.Mat.Structure | None) – Either same, different, or a subset of the non-zero sparsity pattern.
globalup (bool) – Whether global matrices must be updated or not.
- Return type:
Notes
This call makes sense only for spectrum slicing runs, that is, when an interval has been given with
setInterval()andSINVERTis set. And is relevant only when the number of partitions (setKrylovSchurPartitions()) is larger than one.This function modifies the eigenproblem matrices at subcommunicator level, and optionally updates the global matrices in the parent communicator. The updates are expressed as \(A \leftarrow s A + a Au\), \(B \leftarrow t B + b Bu\).
It is possible to update one of the matrices, or both.
The matrices
AuandBumust be equal in all subcommunicators.The
structureflag is passed to thepetsc4py.PETSc.Mat.axpyoperations to perform the updates.If
globalupisTrue, communication is carried out to reconstruct the updated matrices in the parent communicator.
- valuesView(viewer=None)#
Display the computed eigenvalues in a viewer.
Collective.
- Parameters:
viewer (Viewer | None) – Visualization context; if not provided, the standard output is used.
- Return type:
See also
- vectorsView(viewer=None)#
Output computed eigenvectors to a viewer.
Collective.
- Parameters:
viewer (Viewer | None) – Visualization context; if not provided, the standard output is used.
- Return type:
See also
- view(viewer=None)#
Print the EPS data structure.
Collective.
- Parameters:
viewer (Viewer | None) – Visualization context; if not provided, the standard output is used.
- Return type:
See also
Attributes Documentation
- extraction#
The type of extraction technique to be employed.
- max_it#
The maximum iteration count.
- problem_type#
The type of the eigenvalue problem.
- purify#
Eigenvector purification.
- target#
The value of the target.
- tol#
The tolerance.
- track_all#
Compute the residual norm of all approximate eigenpairs.
- true_residual#
Compute the true residual explicitly.
- two_sided#
Two-sided that also computes left eigenvectors.
- which#
The portion of the spectrum to be sought.