Package slepc4py :: Module SLEPc :: Class EPS
[hide private]
[frames] | no frames]

Class EPS


EPS
Nested Classes [hide private]
  Balance
EPS type of balancing used for non-Hermitian problems
  CISSExtraction
EPS CISS extraction technique
  CISSQuadRule
EPS CISS quadrature rule
  Conv
EPS convergence test
  ConvergedReason
EPS convergence reasons
  ErrorType
EPS error type to assess accuracy of computed solutions
  Extraction
EPS extraction technique
  LanczosReorthogType
EPS Lanczos reorthogonalization type
  PowerShiftType
EPS Power shift type.
  ProblemType
EPS problem type
  Stop
EPS stopping test
  Type
EPS type
  Which
EPS desired part of spectrum
Instance Methods [hide private]
a new object with type S, a subtype of T
__new__(S, ...)
 
appendOptionsPrefix(self, prefix)
Appends to the prefix used for searching for all EPS options in the database.
 
cancelMonitor(self)
Clears all monitors for an EPS object.
 
computeError(self, int i, etype=None)
Computes the error (based on the residual norm) associated with the i-th computed eigenpair.
 
create(self, comm=None)
Creates the EPS object.
 
destroy(self)
Destroys the EPS object.
 
errorView(self, etype=None, Viewer viewer=None)
Displays the errors associated with the computed solution (as well as the eigenvalues).
 
getArnoldiDelayed(self)
Gets the type of reorthogonalization used during the Arnoldi iteration.
 
getBV(self)
Obtain the basis vector objects associated to the eigensolver.
 
getBalance(self)
Gets the balancing type used by the EPS object, and the associated parameters.
 
getCISSExtraction(self)
Gets the extraction technique used in the CISS solver.
 
getCISSKSPs(self)
Retrieve the array of linear solver objects associated with the CISS solver.
 
getCISSQuadRule(self)
Gets the quadrature rule used in the CISS solver.
 
getCISSRefinement(self)
Gets the values of various refinement parameters in the CISS solver.
 
getCISSSizes(self)
Gets the values of various size parameters in the CISS solver.
 
getCISSThreshold(self)
Gets the values of various threshold parameters in the CISS solver.
 
getCISSUseST(self)
Gets the flag for using the ST object in the CISS solver.
 
getConverged(self)
Gets the number of converged eigenpairs.
 
getConvergedReason(self)
Gets the reason why the solve() iteration was stopped.
 
getConvergenceTest(self)
Return the method used to compute the error estimate used in the convergence test.
 
getDS(self)
Obtain the direct solver associated to the eigensolver.
 
getDimensions(self)
Gets the number of eigenvalues to compute and the dimension of the subspace.
 
getEigenpair(self, int i, Vec Vr=None, Vec Vi=None)
Gets the i-th solution of the eigenproblem as computed by solve().
 
getEigenvalue(self, int i)
Gets the i-th eigenvalue as computed by solve().
 
getEigenvector(self, int i, Vec Vr, Vec Vi=None)
Gets the i-th eigenvector as computed by solve().
 
getErrorEstimate(self, int i)
Returns the error estimate associated to the i-th computed eigenpair.
 
getExtraction(self)
Gets the extraction type used by the EPS object.
 
getGDBOrth(self)
Returns the orthogonalization used in the search subspace in case of generalized Hermitian problems.
 
getGDBlockSize(self)
Gets the number of vectors to be added to the searching space in every iteration.
 
getGDDoubleExpansion(self)
Gets a flag indicating whether the double expansion variant has been activated or not.
 
getGDInitialSize(self)
Gets the initial size of the searching space.
 
getGDKrylovStart(self)
Gets a flag indicating if the search subspace is started with a Krylov basis.
 
getGDRestart(self)
Gets the number of vectors of the search space after restart and the number of vectors saved from the previous iteration.
 
getInterval(self)
Gets the computational interval for spectrum slicing.
 
getInvariantSubspace(self)
Gets an orthonormal basis of the computed invariant subspace.
 
getIterationNumber(self)
Gets the current iteration number.
 
getJDBOrth(self)
Returns the orthogonalization used in the search subspace in case of generalized Hermitian problems.
 
getJDBlockSize(self)
Gets the number of vectors to be added to the searching space in every iteration.
 
getJDConstCorrectionTol(self)
Returns the flag indicating if the dynamic stopping is being used for solving the correction equation.
 
getJDFix(self)
Gets the threshold for changing the target in the correction equation.
 
getJDInitialSize(self)
Gets the initial size of the searching space.
 
getJDKrylovStart(self)
Gets a flag indicating if the search subspace is started with a Krylov basis.
 
getJDRestart(self)
Gets the number of vectors of the search space after restart and the number of vectors saved from the previous iteration.
 
getKrylovSchurDetectZeros(self)
Gets the flag that enforces zero detection in spectrum slicing.
 
getKrylovSchurDimensions(self)
Gets the dimensions used for each subsolve step in case of doing spectrum slicing for a computational interval.
 
getKrylovSchurInertias(self)
Gets the values of the shifts and their corresponding inertias in case of doing spectrum slicing for a computational interval.
 
getKrylovSchurKSP(self)
Retrieve the linear solver object associated with the internal EPS object in case of doing spectrum slicing for a computational interval.
 
getKrylovSchurLocking(self)
Gets the locking flag used in the Krylov-Schur method.
 
getKrylovSchurPartitions(self)
Gets the number of partitions of the communicator in case of spectrum slicing.
 
getKrylovSchurRestart(self)
Gets the restart parameter used in the Krylov-Schur method.
 
getKrylovSchurSubcommInfo(self)
Gets information related to the case of doing spectrum slicing for a computational interval with multiple communicators.
 
getKrylovSchurSubcommMats(self)
Gets the eigenproblem matrices stored internally in the subcommunicator to which the calling process belongs.
 
getKrylovSchurSubcommPairs(self, int i, Vec V)
Gets the i-th eigenpair stored internally in the multi-communicator to which the calling process belongs.
 
getKrylovSchurSubintervals(self)
Returns the points that delimit the subintervals used in spectrum slicing with several partitions.
 
getLOBPCGBlockSize(self)
Gets the block size used in the LOBPCG method.
 
getLOBPCGLocking(self)
Gets the locking flag used in the LOBPCG method.
 
getLOBPCGRestart(self)
Gets the restart parameter used in the LOBPCG method.
 
getLanczosReorthogType(self)
Gets the type of reorthogonalization used during the Lanczos iteration.
 
getLeftEigenvector(self, int i, Vec Wr, Vec Wi=None)
Gets the i-th left eigenvector as computed by solve().
 
getLyapIIRanks(self)
Return the rank values used for the Lyapunov step.
 
getMonitor(self)
Gets the list of monitor functions.
 
getOperators(self)
Gets the matrices associated with the eigenvalue problem.
 
getOptionsPrefix(self)
Gets the prefix used for searching for all EPS options in the database.
 
getPowerShiftType(self)
Gets the type of shifts used during the power iteration.
 
getProblemType(self)
Gets the problem type from the EPS object.
 
getPurify(self)
Returns the flag indicating whether purification is activated or not.
 
getRG(self)
Obtain the region object associated to the eigensolver.
 
getRQCGReset(self)
Gets the reset parameter used in the RQCG method.
 
getST(self)
Obtain the spectral transformation (ST) object associated to the eigensolver object.
 
getStoppingTest(self)
Gets the stopping function.
 
getTarget(self)
Gets the value of the target.
 
getTolerances(self)
Gets the tolerance and maximum iteration count used by the default EPS convergence tests.
 
getTrackAll(self)
Returns the flag indicating whether all residual norms must be computed or not.
 
getTrueResidual(self)
Returns the flag indicating whether true residual must be computed explicitly or not.
 
getTwoSided(self)
Returns the flag indicating whether a two-sided variant of the algorithm is being used or not.
 
getType(self)
Gets the EPS type of this object.
 
getWhichEigenpairs(self)
Returns which portion of the spectrum is to be sought.
 
isGeneralized(self)
Tells whether the EPS object corresponds to a generalized eigenvalue problem.
 
isHermitian(self)
Tells whether the EPS object corresponds to a Hermitian eigenvalue problem.
 
isPositive(self)
Tells whether the EPS object corresponds to an eigenvalue problem type that requires a positive (semi-) definite matrix B.
 
reset(self)
Resets the EPS object.
 
setArbitrarySelection(self, arbitrary, args=None, kargs=None)
Sets a function to look for eigenvalues according to an arbitrary selection criterion.
 
setArnoldiDelayed(self, delayed)
Activates or deactivates delayed reorthogonalization in the Arnoldi iteration.
 
setBV(self, BV bv)
Associates a basis vectors object to the eigensolver.
 
setBalance(self, balance=None, iterations=None, cutoff=None)
Specifies the balancing technique to be employed by the eigensolver, and some parameters associated to it.
 
setCISSExtraction(self, extraction)
Sets the extraction technique used in the CISS solver.
 
setCISSQuadRule(self, quad)
Sets the quadrature rule used in the CISS solver.
 
setCISSRefinement(self, inner=None, blsize=None)
Sets the values of various refinement parameters in the CISS solver.
 
setCISSSizes(self, ip=None, bs=None, ms=None, npart=None, bsmax=None, realmats=False)
Sets the values of various size parameters in the CISS solver.
 
setCISSThreshold(self, delta=None, spur=None)
Sets the values of various threshold parameters in the CISS solver.
 
setCISSUseST(self, usest)
Sets a flag indicating that the CISS solver will use the ST object for the linear solves.
 
setConvergenceTest(self, conv)
Specifies how to compute the error estimate used in the convergence test.
 
setDS(self, DS ds)
Associates a direct solver object to the eigensolver.
 
setDeflationSpace(self, space)
Add vectors to the basis of the deflation space.
 
setDimensions(self, nev=None, ncv=None, mpd=None)
Sets the number of eigenvalues to compute and the dimension of the subspace.
 
setEigenvalueComparison(self, comparison, args=None, kargs=None)
Specifies the eigenvalue comparison function when setWhichEigenpairs() is set to EPS.Which.USER.
 
setExtraction(self, extraction)
Sets the extraction type used by the EPS object.
 
setFromOptions(self)
Sets EPS options from the options database.
 
setGDBOrth(self, borth)
Selects the orthogonalization that will be used in the search subspace in case of generalized Hermitian problems.
 
setGDBlockSize(self, bs)
Sets the number of vectors to be added to the searching space in every iteration.
 
setGDDoubleExpansion(self, doubleexp)
Activate a variant where 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.
 
setGDInitialSize(self, initialsize)
Sets the initial size of the searching space.
 
setGDKrylovStart(self, krylovstart=True)
Activates or deactivates starting the search subspace with a Krylov basis.
 
setGDRestart(self, minv=None, plusk=None)
Sets the number of vectors of the search space after restart and the number of vectors saved from the previous iteration.
 
setInitialSpace(self, space)
Sets the initial space from which the eigensolver starts to iterate.
 
setInterval(self, inta, intb)
Defines the computational interval for spectrum slicing.
 
setJDBOrth(self, borth)
Selects the orthogonalization that will be used in the search subspace in case of generalized Hermitian problems.
 
setJDBlockSize(self, bs)
Sets the number of vectors to be added to the searching space in every iteration.
 
setJDConstCorrectionTol(self, constant)
Deactivates the dynamic stopping criterion that sets the KSP relative tolerance to 0.5**i, where i is the number of EPS iterations from the last converged value.
 
setJDFix(self, fix)
Sets the threshold for changing the target in the correction equation.
 
setJDInitialSize(self, initialsize)
Sets the initial size of the searching space.
 
setJDKrylovStart(self, krylovstart=True)
Activates or deactivates starting the search subspace with a Krylov basis.
 
setJDRestart(self, minv=None, plusk=None)
Sets the number of vectors of the search space after restart and the number of vectors saved from the previous iteration.
 
setKrylovSchurDetectZeros(self, detect)
Sets a flag to enforce detection of zeros during the factorizations throughout the spectrum slicing computation.
 
setKrylovSchurDimensions(self, nev=None, ncv=None, mpd=None)
Sets the dimensions used for each subsolve step in case of doing spectrum slicing for a computational interval.
 
setKrylovSchurLocking(self, lock)
Choose between locking and non-locking variants of the Krylov-Schur method.
 
setKrylovSchurPartitions(self, npart)
Sets the number of partitions for the case of doing spectrum slicing for a computational interval with the communicator split in several sub-communicators.
 
setKrylovSchurRestart(self, keep)
Sets the restart parameter for the Krylov-Schur method, in particular the proportion of basis vectors that must be kept after restart.
 
setKrylovSchurSubintervals(self, subint)
Sets the subinterval boundaries for spectrum slicing with a computational interval.
 
setLOBPCGBlockSize(self, bs)
Sets the block size of the LOBPCG method.
 
setLOBPCGLocking(self, lock)
Choose between locking and non-locking variants of the LOBPCG method.
 
setLOBPCGRestart(self, restart)
Sets the restart parameter for the LOBPCG method.
 
setLanczosReorthogType(self, reorthog)
Sets the type of reorthogonalization used during the Lanczos iteration.
 
setLeftInitialSpace(self, space)
Sets the left initial space from which the eigensolver starts to iterate.
 
setLyapIIRanks(self, rkc=None, rkl=None)
Set the ranks used in the solution of the Lyapunov equation.
 
setMonitor(self, monitor, args=None, kargs=None)
Appends a monitor function to the list of monitors.
 
setOperators(self, Mat A, Mat B=None)
Sets the matrices associated with the eigenvalue problem.
 
setOptionsPrefix(self, prefix)
Sets the prefix used for searching for all EPS options in the database.
 
setPowerShiftType(self, shift)
Sets the type of shifts used during the power iteration.
 
setProblemType(self, problem_type)
Specifies the type of the eigenvalue problem.
 
setPurify(self, purify=True)
Activate or deactivate eigenvector purification.
 
setRG(self, RG rg)
Associates a region object to the eigensolver.
 
setRQCGReset(self, nrest)
Sets the reset parameter of the RQCG iteration.
 
setST(self, ST st)
Associates a spectral transformation object to the eigensolver.
 
setStoppingTest(self, stopping, args=None, kargs=None)
Sets a function to decide when to stop the outer iteration of the eigensolver.
 
setTarget(self, target)
Sets the value of the target.
 
setTolerances(self, tol=None, max_it=None)
Sets the tolerance and maximum iteration count used by the default EPS convergence tests.
 
setTrackAll(self, trackall)
Specifies if the solver must compute the residual of all approximate eigenpairs or not.
 
setTrueResidual(self, trueres)
Specifies if the solver must compute the true residual explicitly or not.
 
setTwoSided(self, twosided)
Sets the solver to use a two-sided variant so that left eigenvectors are also computed.
 
setType(self, eps_type)
Selects the particular solver to be used in the EPS object.
 
setUp(self)
Sets up all the internal data structures necessary for the execution of the eigensolver.
 
setWhichEigenpairs(self, which)
Specifies which portion of the spectrum is to be sought.
 
solve(self)
Solves the eigensystem.
 
updateKrylovSchurSubcommMats(self, s=1.0, a=1.0, Mat Au=None, t=1.0, b=1.0, Mat Bu=None, structure=None, globalup=False)
Update the eigenproblem matrices stored internally in the subcommunicator to which the calling process belongs.
 
valuesView(self, Viewer viewer=None)
Displays the computed eigenvalues in a viewer.
 
vectorsView(self, Viewer viewer=None)
Outputs computed eigenvectors to a viewer.
 
view(self, Viewer viewer=None)
Prints the EPS data structure.

Inherited from petsc4py.PETSc.Object: __copy__, __deepcopy__, __eq__, __ge__, __gt__, __le__, __lt__, __ne__, __nonzero__, compose, decRef, getAttr, getClassId, getClassName, getComm, getDict, getName, getRefCount, getTabLevel, incRef, incrementTabLevel, query, setAttr, setName, setTabLevel, stateGet, stateIncrease, stateSet, viewFromOptions

Inherited from object: __delattr__, __format__, __getattribute__, __hash__, __init__, __reduce__, __reduce_ex__, __repr__, __setattr__, __sizeof__, __str__, __subclasshook__

Properties [hide private]
  bv
  ds
  extraction
  max_it
  problem_type
  purify
  rg
  st
  target
  tol
  track_all
  true_residual
  two_sided
  which

Inherited from petsc4py.PETSc.Object: classid, comm, fortran, handle, klass, name, prefix, refcount, type

Inherited from object: __class__

Method Details [hide private]

__new__(S, ...)

 
Returns: a new object with type S, a subtype of T
Overrides: object.__new__

appendOptionsPrefix(self, prefix)

 

Appends to the prefix used for searching for all EPS options in the database.

Parameters

prefix: string
The prefix string to prepend to all EPS option requests.
Overrides: petsc4py.PETSc.Object.appendOptionsPrefix

computeError(self, int i, etype=None)

 

Computes the error (based on the residual norm) associated with the i-th computed eigenpair.

Parameters

i: int
Index of the solution to be considered.
etype: EPS.ErrorType enumerate
The error type to compute.

Returns

e: real
The error bound, computed in various ways from the residual norm ||Ax-kBx||_2 where k is the eigenvalue and x is the eigenvector.

Notes

The index i should be a value between 0 and nconv-1 (see getConverged()).

create(self, comm=None)

 

Creates the EPS object.

Parameters

comm: MPI_Comm, optional
MPI communicator; if not provided, it defaults to all processes.

destroy(self)

 
Destroys the EPS object.
Overrides: petsc4py.PETSc.Object.destroy

errorView(self, etype=None, Viewer viewer=None)

 

Displays the errors associated with the computed solution (as well as the eigenvalues).

Parameters

etype: EPS.ErrorType enumerate, optional
The error type to compute.
viewer: Viewer, optional.
Visualization context; if not provided, the standard output is used.

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_DETAIL then a table with eigenvalues and corresponding errors is printed.

getArnoldiDelayed(self)

 

Gets the type of reorthogonalization used during the Arnoldi iteration.

Returns

delayed: bool
True if delayed reorthogonalization is to be used.

getBV(self)

 

Obtain the basis vector objects associated to the eigensolver.

Returns

bv: BV
The basis vectors context.

getBalance(self)

 

Gets the balancing type used by the EPS object, and the associated parameters.

Returns

balance: EPS.Balance enumerate
The balancing method
iterations: int
Number of iterations of the balancing algorithm
cutoff: real
Cutoff value

getCISSExtraction(self)

 

Gets the extraction technique used in the CISS solver.

Returns

extraction: EPS.CISSExtraction enumerate
The extraction technique.

getCISSKSPs(self)

 

Retrieve the array of linear solver objects associated with the CISS solver.

Returns

ksp: list of KSP
The linear solver objects.

Notes

The number of KSP solvers 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.

getCISSQuadRule(self)

 

Gets the quadrature rule used in the CISS solver.

Returns

quad: EPS.CISSQuadRule enumerate
The quadrature rule.

getCISSRefinement(self)

 

Gets the values of various refinement parameters in the CISS solver.

Returns

inner: int
Number of iterative refinement iterations (inner loop).
blsize: int
Number of iterative refinement iterations (blocksize loop).

getCISSSizes(self)

 

Gets the values of various size parameters in the CISS solver.

Returns

ip: int
Number of integration points.
bs: int
Block size.
ms: int
Moment size.
npart: int
Number of partitions when splitting the communicator.
bsmax: int
Maximum block size.
realmats: bool
True if A and B are real.

getCISSThreshold(self)

 

Gets the values of various threshold parameters in the CISS solver.

Returns

delta: float
Threshold for numerical rank.
spur: float
Spurious threshold (to discard spurious eigenpairs.

getCISSUseST(self)

 

Gets the flag for using the ST object in the CISS solver.

Returns

usest: bool
Whether to use the ST object or not.

getConverged(self)

 

Gets the number of converged eigenpairs.

Returns

nconv: int
Number of converged eigenpairs.

Notes

This function should be called after solve() has finished.

getConvergedReason(self)

 

Gets the reason why the solve() iteration was stopped.

Returns

reason: EPS.ConvergedReason enumerate
Negative value indicates diverged, positive value converged.

getConvergenceTest(self)

 

Return the method used to compute the error estimate used in the convergence test.

Returns

conv: EPS.Conv
The method used to compute the error estimate used in the convergence test.

getDS(self)

 

Obtain the direct solver associated to the eigensolver.

Returns

ds: DS
The direct solver context.

getDimensions(self)

 

Gets the number of eigenvalues to compute and the dimension of the subspace.

Returns

nev: int
Number of eigenvalues to compute.
ncv: int
Maximum dimension of the subspace to be used by the solver.
mpd: int
Maximum dimension allowed for the projected problem.

getEigenpair(self, int i, Vec Vr=None, Vec Vi=None)

 

Gets the i-th solution of the eigenproblem as computed by solve(). The solution consists of both the eigenvalue and the eigenvector.

Parameters

i: int
Index of the solution to be obtained.
Vr: Vec
Placeholder for the returned eigenvector (real part).
Vi: Vec
Placeholder for the returned eigenvector (imaginary part).

Returns

e: scalar (possibly complex)
The computed eigenvalue.

Notes

The index i should be a value between 0 and nconv-1 (see getConverged()). Eigenpairs are indexed according to the ordering criterion established with setWhichEigenpairs().

getEigenvalue(self, int i)

 

Gets the i-th eigenvalue as computed by solve().

Parameters

i: int
Index of the solution to be obtained.

Returns

e: scalar (possibly complex)
The computed eigenvalue.

Notes

The index i should be a value between 0 and nconv-1 (see getConverged()). Eigenpairs are indexed according to the ordering criterion established with setWhichEigenpairs().

getEigenvector(self, int i, Vec Vr, Vec Vi=None)

 

Gets the i-th eigenvector as computed by solve().

Parameters

i: int
Index of the solution to be obtained.
Vr: Vec
Placeholder for the returned eigenvector (real part).
Vi: Vec, optional
Placeholder for the returned eigenvector (imaginary part).

Notes

The index i should be a value between 0 and nconv-1 (see getConverged()). Eigenpairs are indexed according to the ordering criterion established with setWhichEigenpairs().

getErrorEstimate(self, int i)

 

Returns the error estimate associated to the i-th computed eigenpair.

Parameters

i: int
Index of the solution to be considered.

Returns

e: real
Error estimate.

Notes

This is the error estimate used internally by the eigensolver. The actual error bound can be computed with computeError().

getExtraction(self)

 

Gets the extraction type used by the EPS object.

Returns

extraction: EPS.Extraction enumerate
The method of extraction.

getGDBOrth(self)

 

Returns the orthogonalization used in the search subspace in case of generalized Hermitian problems.

Returns

borth: bool
Whether to B-orthogonalize the search subspace.

getGDBlockSize(self)

 

Gets the number of vectors to be added to the searching space in every iteration.

Returns

bs: int
The number of vectors added to the search space in every iteration.

getGDDoubleExpansion(self)

 

Gets a flag indicating whether the double expansion variant has been activated or not.

Returns

doubleexp: bool
True if using double expansion.

getGDInitialSize(self)

 

Gets the initial size of the searching space.

Returns

initialsize: int
The number of vectors of the initial searching subspace.

getGDKrylovStart(self)

 

Gets a flag indicating if the search subspace is started with a Krylov basis.

Returns

krylovstart: bool
True if starting the search subspace with a Krylov basis.

getGDRestart(self)

 

Gets the number of vectors of the search space after restart and the number of vectors saved from the previous iteration.

Returns

minv: int
The number of vectors of the search subspace after restart.
plusk: int
The number of vectors saved from the previous iteration.

getInterval(self)

 

Gets the computational interval for spectrum slicing.

Returns

inta: float
The left end of the interval.
intb: float
The right end of the interval.

Notes

If the interval was not set by the user, then zeros are returned.

getInvariantSubspace(self)

 

Gets an orthonormal basis of the computed invariant subspace.

Returns

subspace: list of Vec
Basis of the invariant subspace.

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).

getIterationNumber(self)

 

Gets the current iteration number. If the call to solve() is complete, then it returns the number of iterations carried out by the solution method.

Returns

its: int
Iteration number.

getJDBOrth(self)

 

Returns the orthogonalization used in the search subspace in case of generalized Hermitian problems.

Returns

borth: bool
Whether to B-orthogonalize the search subspace.

getJDBlockSize(self)

 

Gets the number of vectors to be added to the searching space in every iteration.

Returns

bs: int
The number of vectors added to the search space in every iteration.

getJDConstCorrectionTol(self)

 

Returns the flag indicating if the dynamic stopping is being used for solving the correction equation.

Returns

constant: bool
Flag indicating if the dynamic stopping criterion is not being used.

getJDFix(self)

 

Gets the threshold for changing the target in the correction equation.

Returns

fix: float
The threshold for changing the target.

getJDInitialSize(self)

 

Gets the initial size of the searching space.

Returns

initialsize: int
The number of vectors of the initial searching subspace.

getJDKrylovStart(self)

 

Gets a flag indicating if the search subspace is started with a Krylov basis.

Returns

krylovstart: bool
True if starting the search subspace with a Krylov basis.

getJDRestart(self)

 

Gets the number of vectors of the search space after restart and the number of vectors saved from the previous iteration.

Returns

minv: int
The number of vectors of the search subspace after restart.
plusk: int
The number of vectors saved from the previous iteration.

getKrylovSchurDetectZeros(self)

 

Gets the flag that enforces zero detection in spectrum slicing.

Returns

detect: bool
The zero detection flag.

getKrylovSchurDimensions(self)

 

Gets the dimensions used for each subsolve step in case of doing spectrum slicing for a computational interval.

Returns

nev: int
Number of eigenvalues to compute.
ncv: int
Maximum dimension of the subspace to be used by the solver.
mpd: int
Maximum dimension allowed for the projected problem.

getKrylovSchurInertias(self)

 

Gets the values of the shifts and their corresponding inertias in case of doing spectrum slicing for a computational interval.

Returns

shifts: list of float
The values of the shifts used internally in the solver.
inertias: list of int
The values of the inertia in each shift.

getKrylovSchurKSP(self)

 

Retrieve the linear solver object associated with the internal EPS object in case of doing spectrum slicing for a computational interval.

Returns

ksp: KSP
The linear solver object.

getKrylovSchurLocking(self)

 

Gets the locking flag used in the Krylov-Schur method.

Returns

lock: bool
The locking flag.

getKrylovSchurPartitions(self)

 

Gets the number of partitions of the communicator in case of spectrum slicing.

Returns

npart: int
The number of partitions.

getKrylovSchurRestart(self)

 

Gets the restart parameter used in the Krylov-Schur method.

Returns

keep: float
The number of vectors to be kept at restart.

getKrylovSchurSubcommInfo(self)

 

Gets information related to the case of doing spectrum slicing for a computational interval with multiple communicators.

Returns

k: int
Number of the subinterval for the calling process.
n: int
Number of eigenvalues found in the k-th subinterval.
v: Vec
A vector owned by processes in the subcommunicator with dimensions compatible for locally computed eigenvectors.

Notes

This function is only available for spectrum slicing runs.

The returned Vec should be destroyed by the user.

getKrylovSchurSubcommMats(self)

 

Gets the eigenproblem matrices stored internally in the subcommunicator to which the calling process belongs.

Returns

A: Mat
The matrix associated with the eigensystem.
B: Mat
The second matrix in the case of generalized eigenproblems.

Notes

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(self, int i, Vec V)

 

Gets the i-th eigenpair stored internally in the multi-communicator to which the calling process belongs.

Parameters

i: int
Index of the solution to be obtained.
V: Vec
Placeholder for the returned eigenvector.

Returns

e: scalar
The computed eigenvalue.

Notes

The index i should be a value between 0 and n-1, where n is the number of vectors in the local subinterval, see getKrylovSchurSubcommInfo().

getKrylovSchurSubintervals(self)

 

Returns the points that delimit the subintervals used in spectrum slicing with several partitions.

Returns

subint: list of float
Real values specifying subintervals

getLOBPCGBlockSize(self)

 

Gets the block size used in the LOBPCG method.

Returns

bs: int
The block size.

getLOBPCGLocking(self)

 

Gets the locking flag used in the LOBPCG method.

Returns

lock: bool
The locking flag.

getLOBPCGRestart(self)

 

Gets the restart parameter used in the LOBPCG method.

Returns

restart: float
The restart parameter.

getLanczosReorthogType(self)

 

Gets the type of reorthogonalization used during the Lanczos iteration.

Returns

reorthog: EPS.LanczosReorthogType enumerate
The type of reorthogonalization.

getLeftEigenvector(self, int i, Vec Wr, Vec Wi=None)

 

Gets the i-th left eigenvector as computed by solve().

Parameters

i: int
Index of the solution to be obtained.
Wr: Vec
Placeholder for the returned eigenvector (real part).
Wi: Vec, optional
Placeholder for the returned eigenvector (imaginary part).

Notes

The index i should be a value between 0 and nconv-1 (see getConverged()). Eigensolutions are indexed according to the ordering criterion established with setWhichEigenpairs().

Left eigenvectors are available only if the twosided flag was set with setTwoSided().

getLyapIIRanks(self)

 

Return the rank values used for the Lyapunov step.

Returns

rkc: int
The compressed rank.
rkl: int
The Lyapunov rank.

getOperators(self)

 

Gets the matrices associated with the eigenvalue problem.

Returns

A: Mat
The matrix associated with the eigensystem.
B: Mat
The second matrix in the case of generalized eigenproblems.

getOptionsPrefix(self)

 

Gets the prefix used for searching for all EPS options in the database.

Returns

prefix: string
The prefix string set for this EPS object.
Overrides: petsc4py.PETSc.Object.getOptionsPrefix

getPowerShiftType(self)

 

Gets the type of shifts used during the power iteration.

Returns

shift: EPS.PowerShiftType enumerate
The type of shift.

getProblemType(self)

 

Gets the problem type from the EPS object.

Returns

problem_type: EPS.ProblemType enumerate
The problem type that was previously set.

getPurify(self)

 

Returns the flag indicating whether purification is activated or not.

Returns

purify: bool
Whether purification is activated or not.

getRG(self)

 

Obtain the region object associated to the eigensolver.

Returns

rg: RG
The region context.

getRQCGReset(self)

 

Gets the reset parameter used in the RQCG method.

Returns

nrest: int
The number of iterations between resets.

getST(self)

 

Obtain the spectral transformation (ST) object associated to the eigensolver object.

Returns

st: ST
The spectral transformation.

getTarget(self)

 

Gets the value of the target.

Returns

target: float (real or complex)
The value of the target.

Notes

If the target was not set by the user, then zero is returned.

getTolerances(self)

 

Gets the tolerance and maximum iteration count used by the default EPS convergence tests.

Returns

tol: float
The convergence tolerance.
max_it: int
The maximum number of iterations

getTrackAll(self)

 

Returns the flag indicating whether all residual norms must be computed or not.

Returns

trackall: bool
Whether the solver compute all residuals or not.

getTrueResidual(self)

 

Returns the flag indicating whether true residual must be computed explicitly or not.

Returns

trueres: bool
Whether the solver compute all residuals or not.

getTwoSided(self)

 

Returns the flag indicating whether a two-sided variant of the algorithm is being used or not.

Returns

twosided: bool
Whether the two-sided variant is to be used or not.

getType(self)

 

Gets the EPS type of this object.

Returns

type: EPS.Type enumerate
The solver currently being used.
Overrides: petsc4py.PETSc.Object.getType

getWhichEigenpairs(self)

 

Returns which portion of the spectrum is to be sought.

Returns

which: EPS.Which enumerate
The portion of the spectrum to be sought by the solver.

isGeneralized(self)

 

Tells whether the EPS object corresponds to a generalized eigenvalue problem.

Returns

flag: bool
True if two matrices were set with setOperators().

isHermitian(self)

 

Tells whether the EPS object corresponds to a Hermitian eigenvalue problem.

Returns

flag: bool
True if the problem type set with setProblemType() was Hermitian.

isPositive(self)

 

Tells whether the EPS object corresponds to an eigenvalue problem type that requires a positive (semi-) definite matrix B.

Returns

flag: bool
True if the problem type set with setProblemType() was positive.

setArbitrarySelection(self, arbitrary, args=None, kargs=None)

 
Sets 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(self, delayed)

 

Activates or deactivates delayed reorthogonalization in the Arnoldi iteration.

Parameters

delayed: bool
True if delayed reorthogonalization is to be used.

Notes

This call is only relevant if the type was set to EPS.Type.ARNOLDI with setType().

Delayed reorthogonalization is an aggressive optimization for the Arnoldi eigensolver than may provide better scalability, but sometimes makes the solver converge less than the default algorithm.

setBV(self, BV bv)

 

Associates a basis vectors object to the eigensolver.

Parameters

bv: BV
The basis vectors context.

setBalance(self, balance=None, iterations=None, cutoff=None)

 

Specifies the balancing technique to be employed by the eigensolver, and some parameters associated to it.

Parameters

balance: EPS.Balance enumerate
The balancing method
iterations: int
Number of iterations of the balancing algorithm
cutoff: real
Cutoff value

setCISSExtraction(self, extraction)

 

Sets the extraction technique used in the CISS solver.

Parameters

extraction: EPS.CISSExtraction enumerate
The extraction technique.

setCISSQuadRule(self, quad)

 

Sets the quadrature rule used in the CISS solver.

Parameters

quad: EPS.CISSQuadRule enumerate
The quadrature rule.

setCISSRefinement(self, inner=None, blsize=None)

 

Sets the values of various refinement parameters in the CISS solver.

Parameters

inner: int, optional
Number of iterative refinement iterations (inner loop).
blsize: int, optional
Number of iterative refinement iterations (blocksize loop).

setCISSSizes(self, ip=None, bs=None, ms=None, npart=None, bsmax=None, realmats=False)

 

Sets the values of various size parameters in the CISS solver.

Parameters

ip: int, optional
Number of integration points.
bs: int, optional
Block size.
ms: int, optional
Moment size.
npart: int, optional
Number of partitions when splitting the communicator.
bsmax: int, optional
Maximum block size.
realmats: bool, optional
True if A and B are real.

Notes

The default number of partitions is 1. This means the internal KSP object is shared among all processes of the EPS communicator. Otherwise, the communicator is split into npart communicators, so that npart KSP solves proceed simultaneously.

setCISSThreshold(self, delta=None, spur=None)

 

Sets the values of various threshold parameters in the CISS solver.

Parameters

delta: float
Threshold for numerical rank.
spur: float
Spurious threshold (to discard spurious eigenpairs).

setCISSUseST(self, usest)

 

Sets a flag indicating that the CISS solver will use the ST object for the linear solves.

Parameters

usest: bool
Whether to use the ST object or not.

setConvergenceTest(self, conv)

 

Specifies how to compute the error estimate used in the convergence test.

Parameters

conv: EPS.Conv
The method used to compute the error estimate used in the convergence test.

setDS(self, DS ds)

 

Associates a direct solver object to the eigensolver.

Parameters

ds: DS
The direct solver context.

setDeflationSpace(self, space)

 

Add vectors to the basis of the deflation space.

Parameters

space: a Vec or an array of Vec
Set of basis vectors to be added to the deflation space.

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).

The vectors do not need to be mutually orthonormal, since they are explicitly orthonormalized internally.

These vectors do not persist from one solve() call to the other, so the deflation space should be set every time.

setDimensions(self, nev=None, ncv=None, mpd=None)

 

Sets the number of eigenvalues to compute and the dimension of the subspace.

Parameters

nev: int, optional
Number of eigenvalues to compute.
ncv: int, optional
Maximum dimension of the subspace to be used by the solver.
mpd: int, optional
Maximum dimension allowed for the projected problem.

Notes

Use DECIDE for ncv and mpd to assign a reasonably good value, which is dependent on the solution method.

The parameters ncv and mpd are intimately related, so that the user is advised to set one of them at most. Normal usage is the following:

  • In cases where nev is small, the user sets ncv (a reasonable default is 2 * nev).
  • In cases where nev is large, the user sets mpd.

The value of ncv should always be between nev and (nev + mpd), typically ncv = nev + mpd. If nev is not too large, mpd = nev is a reasonable choice, otherwise a smaller value should be used.

setExtraction(self, extraction)

 

Sets the extraction type used by the EPS object.

Parameters

extraction: EPS.Extraction enumerate
The extraction method to be used by the solver.

Notes

Not all eigensolvers support all types of extraction. See the SLEPc documentation for details.

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().

setFromOptions(self)

 

Sets EPS options from the options database. This routine must be called before setUp() if the user is to be allowed to set the solver type.

Notes

To see all options, run your program with the -help option.

Overrides: petsc4py.PETSc.Object.setFromOptions

setGDBOrth(self, borth)

 

Selects the orthogonalization that will be used in the search subspace in case of generalized Hermitian problems.

Parameters

borth: bool
Whether to B-orthogonalize the search subspace.

setGDBlockSize(self, bs)

 

Sets 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.

setGDDoubleExpansion(self, doubleexp)

 

Activate a variant where 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.

Parameters

doubleexp: bool
True if using double expansion.

setGDInitialSize(self, initialsize)

 

Sets the initial size of the searching space.

Parameters

initialsize: int
The number of vectors of the initial searching subspace.

setGDKrylovStart(self, krylovstart=True)

 

Activates or deactivates starting the search subspace with a Krylov basis.

Parameters

krylovstart: bool
True if starting the search subspace with a Krylov basis.

setGDRestart(self, minv=None, plusk=None)

 

Sets the number of vectors of the search space after restart and the number of vectors saved from the previous iteration.

Parameters

minv: int, optional
The number of vectors of the search subspace after restart.
plusk: int, optional
The number of vectors saved from the previous iteration.

setInitialSpace(self, space)

 

Sets the initial space from which the eigensolver starts to iterate.

Parameters

space: Vec or sequence of Vec
The initial space

Notes

Some solvers start to iterate on a single vector (initial vector). In that case, the other vectors are ignored.

In contrast to setDeflationSpace(), 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(self, inta, intb)

 

Defines the computational interval for spectrum slicing.

Parameters

inta: float
The left end of the interval.
intb: float
The right end of the interval.

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, see setWhichEigenpairs().

setJDBOrth(self, borth)

 

Selects the orthogonalization that will be used in the search subspace in case of generalized Hermitian problems.

Parameters

borth: bool
Whether to B-orthogonalize the search subspace.

setJDBlockSize(self, bs)

 

Sets 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.

setJDConstCorrectionTol(self, constant)

 

Deactivates the dynamic stopping criterion that sets the KSP relative tolerance to 0.5**i, where i is the number of EPS iterations from the last converged value.

Parameters

constant: bool
If False, the KSP relative tolerance is set to 0.5**i.

setJDFix(self, fix)

 

Sets the threshold for changing the target in the correction equation.

Parameters

fix: float
The threshold for changing the target.

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.

setJDInitialSize(self, initialsize)

 

Sets the initial size of the searching space.

Parameters

initialsize: int
The number of vectors of the initial searching subspace.

setJDKrylovStart(self, krylovstart=True)

 

Activates or deactivates starting the search subspace with a Krylov basis.

Parameters

krylovstart: bool
True if starting the search subspace with a Krylov basis.

setJDRestart(self, minv=None, plusk=None)

 

Sets the number of vectors of the search space after restart and the number of vectors saved from the previous iteration.

Parameters

minv: int, optional
The number of vectors of the search subspace after restart.
plusk: int, optional
The number of vectors saved from the previous iteration.

setKrylovSchurDetectZeros(self, detect)

 

Sets a flag to enforce detection of zeros during the factorizations throughout the spectrum slicing computation.

Parameters

detect: bool
True if zeros must checked for.

Notes

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(self, nev=None, ncv=None, mpd=None)

 

Sets the dimensions used for each subsolve step in case of doing spectrum slicing for a computational interval. The meaning of the parameters is the same as in setDimensions().

Parameters

nev: int, optional
Number of eigenvalues to compute.
ncv: int, optional
Maximum dimension of the subspace to be used by the solver.
mpd: int, optional
Maximum dimension allowed for the projected problem.

setKrylovSchurLocking(self, lock)

 

Choose between locking and non-locking variants of the Krylov-Schur method.

Parameters

lock: bool
True if the locking variant must be selected.

Notes

The default is to lock converged eigenpairs when the method restarts. This behaviour 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(self, npart)

 

Sets the number of partitions for the case of doing spectrum slicing for a computational interval with the communicator split in several sub-communicators.

Parameters

npart: int
The number of partitions.

Notes

By default, npart=1 so all processes in the communicator participate in the processing of the whole interval. If npart>1 then the interval is divided into npart subintervals, each of them being processed by a subset of processes.

setKrylovSchurRestart(self, keep)

 

Sets the restart parameter for the Krylov-Schur method, in particular the proportion of basis vectors that must be kept after restart.

Parameters

keep: float
The number of vectors to be kept at restart.

Notes

Allowed values are in the range [0.1,0.9]. The default is 0.5.

setKrylovSchurSubintervals(self, subint)

 

Sets the subinterval boundaries for spectrum slicing with a computational interval.

Parameters

subint: list of float
Real values specifying subintervals

Notes

This function must be called after setKrylovSchurPartitions(). For npart partitions, the argument subint must contain npart+1 real values sorted in ascending order: subint_0, subint_1, ..., subint_npart, where the first and last values must coincide with the interval endpoints set with EPSSetInterval(). The subintervals are then defined by two consecutive points: [subint_0,subint_1], [subint_1,subint_2], and so on.

setLOBPCGBlockSize(self, bs)

 

Sets the block size of the LOBPCG method.

Parameters

bs: int
The block size.

setLOBPCGLocking(self, lock)

 

Choose between locking and non-locking variants of the LOBPCG method.

Parameters

lock: bool
True if the locking variant must be selected.

Notes

This flag refers to soft locking (converged vectors within the current block iterate), since hard locking is always used (when nev is larger than the block size).

setLOBPCGRestart(self, restart)

 

Sets the restart parameter for the LOBPCG method. 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.

Parameters

restart: float
The percentage of the block of vectors to force a restart.

Notes

Allowed values are in the range [0.1,1.0]. The default is 0.9.

setLanczosReorthogType(self, reorthog)

 

Sets the type of reorthogonalization used during the Lanczos iteration.

Parameters

reorthog: EPS.LanczosReorthogType enumerate
The type of reorthogonalization.

Notes

This call is only relevant if the type was set to EPS.Type.LANCZOS with setType().

setLeftInitialSpace(self, space)

 

Sets the left initial space from which the eigensolver starts to iterate.

Parameters

space: Vec or sequence of Vec
The left initial space

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.

setLyapIIRanks(self, rkc=None, rkl=None)

 

Set the ranks used in the solution of the Lyapunov equation.

Parameters

rkc: int, optional
The compressed rank.
rkl: int, optional
The Lyapunov rank.

setOperators(self, Mat A, Mat B=None)

 

Sets the matrices associated with the eigenvalue problem.

Parameters

A: Mat
The matrix associated with the eigensystem.
B: Mat, optional
The second matrix in the case of generalized eigenproblems; if not provided, a standard eigenproblem is assumed.

setOptionsPrefix(self, prefix)

 

Sets the prefix used for searching for all EPS options in the database.

Parameters

prefix: string
The prefix string to prepend to all EPS option requests.

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_")
Overrides: petsc4py.PETSc.Object.setOptionsPrefix

setPowerShiftType(self, shift)

 

Sets the type of shifts used during the power iteration. This can be used to emulate the Rayleigh Quotient Iteration (RQI) method.

Parameters

shift: EPS.PowerShiftType enumerate
The type of shift.

Notes

This call is only relevant if the type was set to EPS.Type.POWER with setType().

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.RAYLEIGH or EPS.PowerShiftType.WILKINSON). In this case, the iteration behaves rather like a cubic converging method as RQI.

setProblemType(self, problem_type)

 

Specifies the type of the eigenvalue problem.

Parameters

problem_type: EPS.ProblemType enumerate
The problem type to be set.

Notes

Allowed values are: Hermitian (HEP), non-Hermitian (NHEP), generalized Hermitian (GHEP), generalized non-Hermitian (GNHEP), and generalized non-Hermitian with positive semi-definite B (PGNHEP).

This function must be used to instruct SLEPc to exploit symmetry. 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^H) or generalized Hermitian (i.e. A=A^H, B=B^H, and B positive definite) then it is recommended to set the problem type so that eigensolver can exploit these properties.

setPurify(self, purify=True)

 

Activate or deactivate eigenvector purification.

Parameters

purify: bool, optional
True to activate purification (default).

setRG(self, RG rg)

 

Associates a region object to the eigensolver.

Parameters

rg: RG
The region context.

setRQCGReset(self, nrest)

 

Sets the reset parameter of the RQCG iteration. Every nrest iterations, the solver performs a Rayleigh-Ritz projection step.

Parameters

nrest: int
The number of iterations between resets.

setST(self, ST st)

 

Associates a spectral transformation object to the eigensolver.

Parameters

st: ST
The spectral transformation.

setTarget(self, target)

 

Sets the value of the target.

Parameters

target: float (real or complex)
The value of the target.

Notes

The target is a scalar value used to determine the portion of the spectrum of interest. It is used in combination with setWhichEigenpairs().

setTolerances(self, tol=None, max_it=None)

 

Sets the tolerance and maximum iteration count used by the default EPS convergence tests.

Parameters

tol: float, optional
The convergence tolerance.
max_it: int, optional
The maximum number of iterations

Notes

Use DECIDE for maxits to assign a reasonably good value, which is dependent on the solution method.

setTrackAll(self, trackall)

 

Specifies if the solver must compute the residual of all approximate eigenpairs or not.

Parameters

trackall: bool
Whether compute all residuals or not.

setTrueResidual(self, trueres)

 

Specifies if the solver must compute the true residual explicitly or not.

Parameters

trueres: bool
Whether compute the true residual or not.

setTwoSided(self, twosided)

 

Sets the solver to use a two-sided variant so that left eigenvectors are also computed.

Parameters

twosided: bool
Whether the two-sided variant is to be used or not.

setType(self, eps_type)

 

Selects the particular solver to be used in the EPS object.

Parameters

eps_type: EPS.Type enumerate
The solver to be used.

Notes

See EPS.Type for available methods. The default is EPS.Type.KRYLOVSCHUR. Normally, it is best to use setFromOptions() 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.

setUp(self)

 

Sets up all the internal data structures necessary for the execution of the eigensolver.

Notes

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.

setWhichEigenpairs(self, which)

 

Specifies which portion of the spectrum is to be sought.

Parameters

which: EPS.Which enumerate
The portion of the spectrum to be sought by the solver.

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_IMAGINARY and EPS.Which.SMALLEST_IMAGINARY use the absolute value of the imaginary part for eigenvalue selection.

updateKrylovSchurSubcommMats(self, s=1.0, a=1.0, Mat Au=None, t=1.0, b=1.0, Mat Bu=None, structure=None, globalup=False)

 

Update the eigenproblem matrices stored internally in the subcommunicator to which the calling process belongs.

Parameters

s: float (real or complex)
Scalar that multiplies the existing A matrix.
a: float (real or complex)
Scalar used in the axpy operation on A.
Au: Mat, optional
The matrix used in the axpy operation on A.
t: float (real or complex)
Scalar that multiplies the existing B matrix.
b: float (real or complex)
Scalar used in the axpy operation on B.
Bu: Mat, optional
The matrix used in the axpy operation on B.
structure: PETSc.Mat.Structure enumerate
Either same, different, or a subset of the non-zero sparsity pattern.
globalup: bool
Whether global matrices must be updated or not.

Notes

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 <-- s*A + a*Au, B <-- t*B + b*Bu.

It is possible to update one of the matrices, or both.

The matrices Au and Bu must be equal in all subcommunicators.

The structure flag is passed to the PETSc.Mat.axpy() operations to perform the updates.

If globalup is True, communication is carried out to reconstruct the updated matrices in the parent communicator.

valuesView(self, Viewer viewer=None)

 

Displays the computed eigenvalues in a viewer.

Parameters

viewer: Viewer, optional.
Visualization context; if not provided, the standard output is used.

vectorsView(self, Viewer viewer=None)

 

Outputs computed eigenvectors to a viewer.

Parameters

viewer: Viewer, optional.
Visualization context; if not provided, the standard output is used.

view(self, Viewer viewer=None)

 

Prints the EPS data structure.

Parameters

viewer: Viewer, optional.
Visualization context; if not provided, the standard output is used.
Overrides: petsc4py.PETSc.Object.view