slepc4py.SLEPc.PEP#
- class slepc4py.SLEPc.PEP#
Bases:
ObjectPolynomial Eigenvalue Problem Solver.
The Polynomial Eigenvalue Problem (
PEP) solver is the object provided by slepc4py for specifying a polynomial eigenvalue problem. Apart from the specific solvers for this type of problems, there is anEPS-based solver, i.e., it uses a solver fromEPSto solve a generalized eigenproblem obtained after linearization.Enumerations
PEP basis type for the representation of the polynomial.
PEP CISS extraction technique.
PEP convergence test.
PEP convergence reasons.
PEP error type to assess accuracy of computed solutions.
PEP extraction strategy used.
PEP type of projection to be used in the Jacobi-Davidson solver.
PEP problem type.
PEP refinement strategy.
PEP scheme for solving linear systems during iterative refinement.
PEP scaling strategy.
PEP stopping test.
PEP type.
PEP desired part of spectrum.
Methods Summary
appendOptionsPrefix([prefix])Append to the prefix used for searching for all PEP options in the database.
Clear all monitors for a
PEPobject.computeError(i[, etype])Compute the error associated with the i-th computed eigenpair.
create([comm])Create the PEP object.
destroy()Destroy the PEP object.
errorView([etype, viewer])Display the errors associated with the computed solution.
getBV()Get the basis vectors object associated to the eigensolver.
getBasis()Get the type of polynomial basis used.
Get the extraction technique used in the CISS solver.
Get the array of linear solver objects associated with 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 number of converged eigenpairs.
Get the reason why the
solve()iteration was stopped.Get the method used to compute the error estimate used in the convergence test.
getDS()Get the direct solver associated to the eigensolver.
Get the 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 eigenvalue comparison function.
Get the error estimate associated to the i-th computed eigenpair.
Get the extraction technique used by the
PEPobject.Get the computational interval for spectrum slicing.
Get the current iteration number.
getJDFix()Get threshold for changing the target in the correction equation.
Get the maximum allowed value of the minimality index.
Get the type of projection to be used in the Jacobi-Davidson solver.
Get the restart parameter used in the Jacobi-Davidson method.
Get the flag for reusing the preconditioner.
Get the eigensolver object associated to the polynomial eigenvalue solver.
Get if the matrices for the linearization are built explicitly.
Get the coeffs.
Get the list of monitor functions.
Get the matrices associated with the eigenvalue problem.
Get the prefix used for searching for all PEP options in the database.
Get the problem type from the PEP object.
Get the locking flag used in the Q-Arnoldi method.
Get the restart parameter used in the Q-Arnoldi method.
getRG()Get the region object associated to the eigensolver.
Get the refinement strategy used by the PEP object.
Get the
KSPobject used by the eigensolver in the refinement phase.getST()Get the spectral transformation object associated to the eigensolver.
Get the flag for the eigenvalue type check in spectrum slicing.
Get the flag that enforces zero detection in spectrum slicing.
Get the dimensions used for each subsolve step.
Get the values of the shifts and their corresponding inertias.
Get the coefficients that define the linearization of a quadratic eigenproblem.
Get the locking flag used in the STOAR method.
getScale([Dl, Dr])Get the strategy used for scaling the polynomial eigenproblem.
Get the stopping test function.
Get the locking flag used in the TOAR method.
Get the restart parameter used in the TOAR method.
Get the value of the target.
Get the tolerance and maximum iteration count.
Get the flag indicating whether all residual norms must be computed.
getType()Get the PEP type of this object.
Get which portion of the spectrum is to be sought.
reset()Reset the PEP object.
setBV(bv)Set a basis vectors object associated to the eigensolver.
setBasis(basis)Set the type of polynomial basis used.
setCISSExtraction(extraction)Set the extraction technique 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.
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.
setDimensions([nev, ncv, mpd])Set the number of eigenvalues to compute and the dimension of the subspace.
setEigenvalueComparison(comparison[, args, ...])Set an eigenvalue comparison function.
setExtract(extract)Set the extraction strategy to be used.
Set PEP options from the options database.
setInitialSpace(space)Set the initial space from which the eigensolver starts to iterate.
setInterval(inta, intb)Set the computational interval for spectrum slicing.
setJDFix(fix)Set the threshold for changing the target in the correction equation.
setJDMinimalityIndex(flag)Set the maximum allowed value for the minimality index.
setJDProjection(proj)Set the type of projection to be used in the Jacobi-Davidson solver.
setJDRestart(keep)Set the restart parameter for the Jacobi-Davidson method.
setJDReusePreconditioner(flag)Set a flag indicating whether the preconditioner must be reused or not.
setLinearEPS(eps)Set an eigensolver object associated to the polynomial eigenvalue solver.
setLinearExplicitMatrix(flag)Set flag to explicitly build the matrices for the linearization.
setLinearLinearization([alpha, beta])Set the coefficients that define the linearization of a quadratic eigenproblem.
setMonitor(monitor[, args, kargs])Append a monitor function to the list of monitors.
setOperators(operators)Set the matrices associated with the eigenvalue problem.
setOptionsPrefix([prefix])Set the prefix used for searching for all PEP options in the database.
setProblemType(problem_type)Set the type of the polynomial eigenvalue problem.
setQArnoldiLocking(lock)Toggle between locking and non-locking variants of the Q-Arnoldi method.
setQArnoldiRestart(keep)Set the restart parameter for the Q-Arnoldi method.
setRG(rg)Set a region object associated to the eigensolver.
setRefine(ref[, npart, tol, its, scheme])Set the refinement strategy used by the PEP object.
setST(st)Set a spectral transformation object associated to the eigensolver.
Set flag to check if all eigenvalues have the same definite type.
setSTOARDetectZeros(detect)Set flag to enforce detection of zeros during the factorizations.
setSTOARDimensions([nev, ncv, mpd])Set the dimensions used for each subsolve step.
setSTOARLinearization([alpha, beta])Set the coefficients that define the linearization of a quadratic eigenproblem.
setSTOARLocking(lock)Toggle between locking and non-locking variants of the STOAR method.
setScale(scale[, alpha, Dl, Dr, its, lbda])Set the scaling strategy to be used.
setStoppingTest(stopping[, args, kargs])Set a function to decide when to stop the outer iteration of the eigensolver.
setTOARLocking(lock)Toggle between locking and non-locking variants of the TOAR method.
setTOARRestart(keep)Set the restart parameter for the TOAR method.
setTarget(target)Set the value of the target.
setTolerances([tol, max_it])Set the tolerance and maximum iteration count.
setTrackAll(trackall)Set flag to compute the residual of all approximate eigenpairs.
setType(pep_type)Set the particular solver to be used in the PEP object.
setUp()Set up all the internal data structures.
setWhichEigenpairs(which)Set which portion of the spectrum is to be sought.
solve()Solve the polynomial eigenproblem.
valuesView([viewer])Display the computed eigenvalues in a viewer.
vectorsView([viewer])Output computed eigenvectors to a viewer.
view([viewer])Print the PEP 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.
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.
The portion of the spectrum to be sought.
Methods Documentation
- appendOptionsPrefix(prefix=None)#
Append to the prefix used for searching for all PEP options in the database.
Logically collective.
- cancelMonitor()#
Clear all monitors for a
PEPobject.Logically collective.
See also
Source code at slepc4py/SLEPc/PEP.pyx:1589
- 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 \(\|P(\lambda)x\|_2\) where \(\lambda\) is the eigenvalue and \(x\) is the eigenvector.
- Return type:
Notes
The index
ishould be a value between0andnconv-1(seegetConverged()).See also
- create(comm=None)#
Create the PEP object.
Collective.
- Parameters:
comm (Comm | None) – MPI communicator. If not provided, it defaults to all processes.
- Return type:
See also
- destroy()#
Destroy the PEP object.
Collective.
See also
Source code at slepc4py/SLEPc/PEP.pyx:324
- 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
- getBV()#
Get the basis vectors object associated to the eigensolver.
Not collective.
- Returns:
The basis vectors context.
- Return type:
- getBasis()#
Get the type of polynomial basis used.
Not collective.
- Returns:
The basis that was previously set.
- 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.
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
- 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
- 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 the method used 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 the 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:
The computed eigenvalue.
- Return type:
Notes
The index
ishould be a value between0andnconv-1(seegetConverged()). Eigenpairs are indexed according to the ordering criterion established withsetWhichEigenpairs().The eigenvector is normalized to have unit norm.
See also
- getEigenvalueComparison()#
Get the eigenvalue comparison function.
Not collective.
- Returns:
The eigenvalue comparison function.
- Return type:
See also
- 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
- getExtract()#
Get the extraction technique used by the
PEPobject.Not collective.
- Returns:
The extraction strategy.
- 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
- 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
- getJDFix()#
Get threshold for changing the target in the correction equation.
Not collective.
- Returns:
The threshold for changing the target.
- Return type:
See also
- getJDMinimalityIndex()#
Get the maximum allowed value of the minimality index.
Not collective.
- Returns:
The maximum minimality index.
- Return type:
See also
- getJDProjection()#
Get the type of projection to be used in the Jacobi-Davidson solver.
Not collective.
- Returns:
The type of projection.
- Return type:
See also
- getJDRestart()#
Get the restart parameter used in the Jacobi-Davidson method.
Not collective.
- Returns:
The number of vectors to be kept at restart.
- Return type:
See also
- getJDReusePreconditioner()#
Get the flag for reusing the preconditioner.
Not collective.
- Returns:
The reuse flag.
- Return type:
- getLinearEPS()#
Get the eigensolver object associated to the polynomial eigenvalue solver.
Collective.
- Returns:
The linear eigensolver.
- Return type:
See also
- getLinearExplicitMatrix()#
Get if the matrices for the linearization are built explicitly.
Not collective.
- Returns:
Boolean flag indicating if the matrices are built explicitly.
- Return type:
- getLinearLinearization()#
Get the coeffs. defining the linearization of a quadratic eigenproblem.
Not collective.
- Returns:
- Return type:
- 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:
The matrices associated with the eigensystem.
- Return type:
See also
- getOptionsPrefix()#
Get the prefix used for searching for all PEP options in the database.
Not collective.
- Returns:
The prefix string set for this PEP object.
- Return type:
- getProblemType()#
Get the problem type from the PEP object.
Not collective.
- Returns:
The problem type that was previously set.
- Return type:
See also
- getQArnoldiLocking()#
Get the locking flag used in the Q-Arnoldi method.
Not collective.
- Returns:
The locking flag.
- Return type:
See also
- getQArnoldiRestart()#
Get the restart parameter used in the Q-Arnoldi method.
Not collective.
- Returns:
The number of vectors to be kept at restart.
- Return type:
See also
- getRG()#
Get the region object associated to the eigensolver.
Not collective.
- Returns:
The region context.
- Return type:
- getRefine()#
Get the refinement strategy used by the PEP object.
Not collective.
- Returns:
ref (
Refine) – The refinement type.npart (
int) – The number of partitions of the communicator.tol (
float) – The convergence tolerance.its (
int) – The maximum number of refinement iterations.scheme (
RefineScheme) – Scheme for solving linear systems.
- Return type:
See also
- getRefineKSP()#
Get the
KSPobject used by the eigensolver in the refinement phase.Collective.
- Returns:
The linear solver object.
- Return type:
See also
- getST()#
Get the spectral transformation object associated to the eigensolver.
Not collective.
- Returns:
The spectral transformation.
- Return type:
- getSTOARCheckEigenvalueType()#
Get the flag for the eigenvalue type check in spectrum slicing.
Not collective.
- Returns:
Whether the eigenvalue type is checked or not.
- Return type:
- getSTOARDetectZeros()#
Get the flag that enforces zero detection in spectrum slicing.
Not collective.
- Returns:
The zero detection flag.
- Return type:
See also
- getSTOARDimensions()#
Get the dimensions used for each subsolve step.
Not collective.
- Returns:
- Return type:
See also
- getSTOARInertias()#
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.See also
- getSTOARLinearization()#
Get the coefficients that define the linearization of a quadratic eigenproblem.
Not collective.
- Returns:
- Return type:
- getSTOARLocking()#
Get the locking flag used in the STOAR method.
Not collective.
- Returns:
The locking flag.
- Return type:
See also
- getScale(Dl=None, Dr=None)#
Get the strategy used for scaling the polynomial eigenproblem.
Not collective.
- Parameters:
Dl (petsc4py.PETSc.Vec | None) – Placeholder for the returned left diagonal matrix.
Dr (petsc4py.PETSc.Vec | None) – Placeholder for the returned right diagonal matrix.
- Returns:
- Return type:
See also
- getStoppingTest()#
Get the stopping test function.
Not collective.
- Returns:
The stopping test function.
- Return type:
See also
- getTOARLocking()#
Get the locking flag used in the TOAR method.
Not collective.
- Returns:
The locking flag.
- Return type:
See also
- getTOARRestart()#
Get the restart parameter used in the TOAR method.
Not collective.
- Returns:
The number of vectors to be kept at restart.
- 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
- getTolerances()#
Get the tolerance and maximum iteration count.
Not collective.
Get the tolerance and maximum iteration count used by the default PEP convergence tests.
- Returns:
- Return type:
See also
- getTrackAll()#
Get the flag indicating whether all residual norms must be computed.
Not collective.
- Returns:
Whether the solver computes all residuals or not.
- Return type:
See also
- getType()#
Get the PEP 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
- reset()#
Reset the PEP object.
Collective.
See also
Source code at slepc4py/SLEPc/PEP.pyx:338
- Return type:
- setBV(bv)#
Set a basis vectors object associated to the eigensolver.
Collective.
- setBasis(basis)#
Set the type of polynomial basis used.
Logically collective.
Set the type of polynomial basis used to describe the polynomial eigenvalue problem.
Notes
By default, the coefficient matrices passed via
setOperators()are expressed in the monomial basis, i.e. \(P(\lambda)=A_0+\lambda A_1+\lambda^2 A_2+\dots+\lambda^d A_d\). Other polynomial bases may have better numerical behavior, but the user must then pass the coefficient matrices accordingly.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
- 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 thePEPcommunicator. Otherwise, the communicator is split intonpartcommunicators, so thatnpartpetsc4py.PETSc.KSPsolves proceed simultaneously.See also
getCISSSizes,setCISSThreshold,setCISSRefinement,PEPCISSSetSizes
- setCISSThreshold(delta=None, spur=None)#
Set the values of various threshold parameters in the CISS solver.
Logically collective.
- Parameters:
- Return type:
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.
- setDimensions(nev=None, ncv=None, mpd=None)#
Set the 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.See also
- setEigenvalueComparison(comparison, args=None, kargs=None)#
Set an eigenvalue comparison function.
Logically collective.
Notes
This eigenvalue comparison function is used when
setWhichEigenpairs()is set toPEP.Which.USER.
- setExtract(extract)#
Set the extraction strategy to be used.
Logically collective.
Notes
This is relevant for solvers based on linearization. Once the solver has converged, the polynomial eigenvectors can be extracted from the eigenvectors of the linearized problem in different ways.
See also
- setFromOptions()#
Set PEP 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/PEP.pyx:486
- Return type:
- 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.
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.
See also
- 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 quadratic eigenproblems in a given interval. This function provides the interval to be considered. It must be used in combination with
PEP.Which.ALL, seesetWhichEigenpairs(). Note that in polynomial eigenproblems spectrum slicing is implemented inSTOARonly.See also
- 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 fix value, the target is set to the corresponding eigenvalue.
See also
- setJDMinimalityIndex(flag)#
Set the maximum allowed value for the minimality index.
Logically collective.
Notes
The default value is equal to the degree of the polynomial. A smaller value can be used if the wanted eigenvectors are known to be linearly independent.
See also
- setJDProjection(proj)#
Set the type of projection to be used in the Jacobi-Davidson solver.
Logically collective.
- Parameters:
proj (JDProjection) – The type of projection.
- Return type:
See also
- setJDRestart(keep)#
Set the restart parameter for the Jacobi-Davidson method.
Logically collective.
Set the restart parameter for the Jacobi-Davidson method, in particular 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.
See also
- setJDReusePreconditioner(flag)#
Set a flag indicating whether the preconditioner must be reused or not.
Logically collective.
Notes
The default value is
False. If set toTrue, the preconditioner is built only at the beginning, using the target value. Otherwise, it may be rebuilt (depending on thefixparameter) at each iteration from the Ritz value.
- setLinearEPS(eps)#
Set an eigensolver object associated to the polynomial eigenvalue solver.
Collective.
See also
- setLinearExplicitMatrix(flag)#
Set flag to explicitly build the matrices for the linearization.
Logically collective.
- Parameters:
flag (bool) – Boolean flag indicating if the matrices are built explicitly.
- Return type:
- setLinearLinearization(alpha=1.0, beta=0.0)#
Set the coefficients that define the linearization of a quadratic eigenproblem.
Logically collective.
- Parameters:
- Return type:
- setMonitor(monitor, args=None, kargs=None)#
Append a monitor function to the list of monitors.
Logically collective.
See also
- setOperators(operators)#
Set the matrices associated with the eigenvalue problem.
Collective.
Notes
The polynomial eigenproblem is defined as \(P(\lambda)x=0\), where \(\lambda\) is the eigenvalue, \(x\) is the eigenvector, and \(P\) is defined as \(P(\lambda) = A_0 + \lambda A_1 + \dots + \lambda^d A_d\), with \(d\) =
nmat-1 (the degree of \(P\)). For non-monomial bases, this expression is different.See also
- setOptionsPrefix(prefix=None)#
Set the prefix used for searching for all PEP options in the database.
Logically collective.
- Parameters:
prefix (str | None) – The prefix string to prepend to all PEP 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 PEP contexts, one could call:
P1.setOptionsPrefix("pep1_") P2.setOptionsPrefix("pep2_")
- setProblemType(problem_type)#
Set the type of the polynomial eigenvalue problem.
Logically collective.
- Parameters:
problem_type (ProblemType) – The problem type to be set.
- Return type:
Notes
This function is used to instruct SLEPc to exploit certain structure in the polynomial eigenproblem. By default, no particular structure is assumed.
If the problem matrices are Hermitian (symmetric in the real case) or Hermitian/skew-Hermitian then the solver can exploit this fact to perform less operations or provide better stability. Hyperbolic problems are a particular case of Hermitian problems, some solvers may treat them simply as Hermitian.
See also
- setQArnoldiLocking(lock)#
Toggle between locking and non-locking variants of the Q-Arnoldi 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).
See also
- setQArnoldiRestart(keep)#
Set the restart parameter for the Q-Arnoldi method.
Logically collective.
Set the restart parameter for the Q-Arnoldi method, in particular 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.
See also
- setRG(rg)#
Set a region object associated to the eigensolver.
Collective.
- setRefine(ref, npart=None, tol=None, its=None, scheme=None)#
Set the refinement strategy used by the PEP object.
Logically collective.
Set the refinement strategy used by the PEP object, and the associated parameters.
- Parameters:
- Return type:
See also
- setST(st)#
Set a spectral transformation object associated to the eigensolver.
Collective.
- setSTOARCheckEigenvalueType(flag)#
Set flag to check if all eigenvalues have the same definite type.
Logically collective.
Set a flag to check that all the eigenvalues obtained throughout the spectrum slicing computation have the same definite type.
Notes
This option is relevant only for spectrum slicing computations, but is ignored in
slepc4py.SLEPc.PEP.ProblemType.HYPERBOLICproblems.This flag is turned on by default, to guarantee that the computed eigenvalues have the same type (otherwise the computed solution might be wrong). But since the check is computationally quite expensive, the check may be turned off if the user knows for sure that all eigenvalues in the requested interval have the same type.
- setSTOARDetectZeros(detect)#
Set flag to enforce detection of zeros during the factorizations.
Logically collective.
Set a flag to enforce 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. This feature currently requires an external package for factorizations with support for zero detection, e.g. MUMPS.
See also
- setSTOARDimensions(nev=None, ncv=None, mpd=None)#
Set the dimensions used for each subsolve step.
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 childPEPobject.
- setSTOARLinearization(alpha=1.0, beta=0.0)#
Set the coefficients that define the linearization of a quadratic eigenproblem.
Logically collective.
- Parameters:
- Return type:
- setSTOARLocking(lock)#
Toggle between locking and non-locking variants of the STOAR 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).
See also
- setScale(scale, alpha=None, Dl=None, Dr=None, its=None, lbda=None)#
Set the scaling strategy to be used.
Collective.
- Parameters:
scale (Scale) – The scaling strategy.
Dl (petsc4py.PETSc.Vec | None) – The left diagonal matrix.
Dr (petsc4py.PETSc.Vec | None) – The right diagonal matrix.
its (int | None) – The number of iterations of diagonal scaling.
lbda (float | None) – Approximation of the wanted eigenvalues (modulus).
- Return type:
See also
- setStoppingTest(stopping, args=None, kargs=None)#
Set a function to decide when to stop the outer iteration of the eigensolver.
Logically collective.
See also
- setTOARLocking(lock)#
Toggle between locking and non-locking variants of the TOAR 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).
See also
- setTOARRestart(keep)#
Set the restart parameter for the TOAR method.
Logically collective.
Set the restart parameter for the TOAR method, in particular 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.
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
- setTolerances(tol=None, max_it=None)#
Set the tolerance and maximum iteration count.
Logically collective.
Set the tolerance and maximum iteration count used by the default PEP convergence tests.
- Parameters:
- Return type:
Notes
Use
DETERMINEformax_itto assign a reasonably good value, which is dependent on the solution method.See also
- setTrackAll(trackall)#
Set flag to compute the residual of all approximate eigenpairs.
Logically collective.
Set if the solver must compute the residual of all approximate eigenpairs or not.
See also
- setType(pep_type)#
Set the particular solver to be used in the PEP object.
Logically collective.
Notes
The default is
TOAR. Normally, it is best to usesetFromOptions()and then set the PEP 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.Source code at slepc4py/SLEPc/PEP.pyx:1604
- 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 PEP account for all the possible values. Also, some values make sense only for certain types of problems. If SLEPc is compiled for real numbers
PEP.Which.LARGEST_IMAGINARYandPEP.Which.SMALLEST_IMAGINARYuse the absolute value of the imaginary part for eigenvalue selection.The target is a scalar value provided with
setTarget().The criterion
PEP.Which.TARGET_IMAGINARYis available only in case PETSc and SLEPc have been built with complex scalars.PEP.Which.ALLis intended for use in combination with an interval (seesetInterval()), when all eigenvalues within the interval are requested, or in the context of thePEP.Type.CISSsolver for computing all eigenvalues in a region.See also
getWhichEigenpairs,setTarget,setInterval,PEPSetWhichEigenpairs
- solve()#
Solve the polynomial eigenproblem.
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,PEPSolveSource code at slepc4py/SLEPc/PEP.pyx:1626
- Return type:
- 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 PEP data structure.
Collective.
- Parameters:
viewer (Viewer | None) – Visualization context; if not provided, the standard output is used.
- Return type:
See also
Attributes Documentation
- extract#
The type of extraction technique to be employed.
- max_it#
The maximum iteration count.
- problem_type#
The type of the eigenvalue problem.
- target#
The value of the target.
- tol#
The tolerance.
- track_all#
Compute the residual norm of all approximate eigenpairs.
- which#
The portion of the spectrum to be sought.