slepc4py.SLEPc.PEP#

class slepc4py.SLEPc.PEP#

Bases: Object

Polynomial 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 an EPS-based solver, i.e., it uses a solver from EPS to solve a generalized eigenproblem obtained after linearization.

Enumerations

Basis

PEP basis type for the representation of the polynomial.

CISSExtraction

PEP CISS extraction technique.

Conv

PEP convergence test.

ConvergedReason

PEP convergence reasons.

ErrorType

PEP error type to assess accuracy of computed solutions.

Extract

PEP extraction strategy used.

JDProjection

PEP type of projection to be used in the Jacobi-Davidson solver.

ProblemType

PEP problem type.

Refine

PEP refinement strategy.

RefineScheme

PEP scheme for solving linear systems during iterative refinement.

Scale

PEP scaling strategy.

Stop

PEP stopping test.

Type

PEP type.

Which

PEP desired part of spectrum.

Methods Summary

appendOptionsPrefix([prefix])

Append to the prefix used for searching for all PEP options in the database.

cancelMonitor()

Clear all monitors for a PEP object.

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.

getCISSExtraction()

Get the extraction technique used in the CISS solver.

getCISSKSPs()

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

getCISSRefinement()

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

getCISSSizes()

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

getCISSThreshold()

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

getConverged()

Get the number of converged eigenpairs.

getConvergedReason()

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

getConvergenceTest()

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

getDS()

Get the direct solver associated to the eigensolver.

getDimensions()

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

getEigenvalueComparison()

Get the eigenvalue comparison function.

getErrorEstimate(i)

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

getExtract()

Get the extraction technique used by the PEP object.

getInterval()

Get the computational interval for spectrum slicing.

getIterationNumber()

Get the current iteration number.

getJDFix()

Get threshold for changing the target in the correction equation.

getJDMinimalityIndex()

Get the maximum allowed value of the minimality index.

getJDProjection()

Get the type of projection to be used in the Jacobi-Davidson solver.

getJDRestart()

Get the restart parameter used in the Jacobi-Davidson method.

getJDReusePreconditioner()

Get the flag for reusing the preconditioner.

getLinearEPS()

Get the eigensolver object associated to the polynomial eigenvalue solver.

getLinearExplicitMatrix()

Get if the matrices for the linearization are built explicitly.

getLinearLinearization()

Get the coeffs.

getMonitor()

Get the list of monitor functions.

getOperators()

Get the matrices associated with the eigenvalue problem.

getOptionsPrefix()

Get the prefix used for searching for all PEP options in the database.

getProblemType()

Get the problem type from the PEP object.

getQArnoldiLocking()

Get the locking flag used in the Q-Arnoldi method.

getQArnoldiRestart()

Get the restart parameter used in the Q-Arnoldi method.

getRG()

Get the region object associated to the eigensolver.

getRefine()

Get the refinement strategy used by the PEP object.

getRefineKSP()

Get the KSP object used by the eigensolver in the refinement phase.

getST()

Get the spectral transformation object associated to the eigensolver.

getSTOARCheckEigenvalueType()

Get the flag for the eigenvalue type check in spectrum slicing.

getSTOARDetectZeros()

Get the flag that enforces zero detection in spectrum slicing.

getSTOARDimensions()

Get the dimensions used for each subsolve step.

getSTOARInertias()

Get the values of the shifts and their corresponding inertias.

getSTOARLinearization()

Get the coefficients that define the linearization of a quadratic eigenproblem.

getSTOARLocking()

Get the locking flag used in the STOAR method.

getScale([Dl, Dr])

Get the strategy used for scaling the polynomial eigenproblem.

getStoppingTest()

Get the stopping test function.

getTOARLocking()

Get the locking flag used in the TOAR method.

getTOARRestart()

Get the restart parameter used in the TOAR method.

getTarget()

Get the value of the target.

getTolerances()

Get the tolerance and maximum iteration count.

getTrackAll()

Get the flag indicating whether all residual norms must be computed.

getType()

Get the PEP type of this object.

getWhichEigenpairs()

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.

setFromOptions()

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.

setSTOARCheckEigenvalueType(flag)

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

bv

The basis vectors (BV) object associated.

ds

The direct solver (DS) object associated.

extract

The type of extraction technique to be employed.

max_it

The maximum iteration count.

problem_type

The type of the eigenvalue problem.

rg

The region (RG) object associated.

st

The spectral transformation (ST) object associated.

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.

Methods Documentation

appendOptionsPrefix(prefix=None)#

Append to 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:

None

Source code at slepc4py/SLEPc/PEP.pyx:467

cancelMonitor()#

Clear all monitors for a PEP object.

Logically collective.

See also

PEPMonitorCancel

Source code at slepc4py/SLEPc/PEP.pyx:1589

Return type:

None

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:
  • i (int) – Index of the solution to be considered.

  • etype (ErrorType | None) – The error type to compute.

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:

float

Notes

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

Source code at slepc4py/SLEPc/PEP.pyx:1785

create(comm=None)#

Create the PEP object.

Collective.

Parameters:

comm (Comm | None) – MPI communicator. If not provided, it defaults to all processes.

Return type:

Self

See also

PEPCreate

Source code at slepc4py/SLEPc/PEP.pyx:350

destroy()#

Destroy the PEP object.

Collective.

See also

PEPDestroy

Source code at slepc4py/SLEPc/PEP.pyx:324

Return type:

Self

errorView(etype=None, viewer=None)#

Display the errors associated with the computed solution.

Collective.

Display the errors and the eigenvalues.

Parameters:
Return type:

None

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.

Source code at slepc4py/SLEPc/PEP.pyx:1823

getBV()#

Get the basis vectors object associated to the eigensolver.

Not collective.

Returns:

The basis vectors context.

Return type:

BV

See also

setBV, PEPGetBV

Source code at slepc4py/SLEPc/PEP.pyx:1250

getBasis()#

Get the type of polynomial basis used.

Not collective.

Returns:

The basis that was previously set.

Return type:

Basis

See also

setBasis, PEPGetBasis

Source code at slepc4py/SLEPc/PEP.pyx:505

getCISSExtraction()#

Get the extraction technique used in the CISS solver.

Not collective.

Returns:

The extraction technique.

Return type:

CISSExtraction

Source code at slepc4py/SLEPc/PEP.pyx:2729

getCISSKSPs()#

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

Collective.

Returns:

The linear solver objects.

Return type:

list of petsc4py.PETSc.KSP

Notes

The number of petsc4py.PETSc.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.

Source code at slepc4py/SLEPc/PEP.pyx:2926

getCISSRefinement()#

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

Not collective.

Returns:

  • inner (int) – Number of iterative refinement iterations (inner loop).

  • blsize (int) – Number of iterative refinement iterations (blocksize loop).

Return type:

tuple[int, int]

Source code at slepc4py/SLEPc/PEP.pyx:2904

getCISSSizes()#

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

Not collective.

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.

Return type:

tuple[int, int, int, int, int, bool]

Source code at slepc4py/SLEPc/PEP.pyx:2802

getCISSThreshold()#

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

Not collective.

Returns:

  • delta (float) – Threshold for numerical rank.

  • spur (float) – Spurious threshold (to discard spurious eigenpairs.

Return type:

tuple[float, float]

Source code at slepc4py/SLEPc/PEP.pyx:2859

getConverged()#

Get the number of converged eigenpairs.

Not collective.

Returns:

nconv – Number of converged eigenpairs.

Return type:

int

Notes

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

The value nconv may be different from the number of requested solutions nev, but not larger than ncv, see setDimensions().

Source code at slepc4py/SLEPc/PEP.pyx:1689

getConvergedReason()#

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

Not collective.

Returns:

Negative value indicates diverged, positive value converged.

Return type:

ConvergedReason

Source code at slepc4py/SLEPc/PEP.pyx:1670

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:

Conv

Source code at slepc4py/SLEPc/PEP.pyx:820

getDS()#

Get the direct solver associated to the eigensolver.

Not collective.

Returns:

The direct solver context.

Return type:

DS

See also

setDS, PEPGetDS

Source code at slepc4py/SLEPc/PEP.pyx:1324

getDimensions()#

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

Not collective.

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.

Return type:

tuple[int, int, int]

Source code at slepc4py/SLEPc/PEP.pyx:1038

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:
  • i (int) – Index of the solution to be obtained.

  • Vr (Vec | None) – Placeholder for the returned eigenvector (real part).

  • Vi (Vec | None) – Placeholder for the returned eigenvector (imaginary part).

Returns:

The computed eigenvalue.

Return type:

complex

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

The eigenvector is normalized to have unit norm.

Source code at slepc4py/SLEPc/PEP.pyx:1715

getEigenvalueComparison()#

Get the eigenvalue comparison function.

Not collective.

Returns:

The eigenvalue comparison function.

Return type:

PEPEigenvalueComparison

Source code at slepc4py/SLEPc/PEP.pyx:1530

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:

float

Notes

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

Source code at slepc4py/SLEPc/PEP.pyx:1756

getExtract()#

Get the extraction technique used by the PEP object.

Not collective.

Returns:

The extraction strategy.

Return type:

Extract

Source code at slepc4py/SLEPc/PEP.pyx:979

getInterval()#

Get the computational interval for spectrum slicing.

Not collective.

Returns:

  • inta (float) – The left end of the interval.

  • intb (float) – The right end of the interval.

Return type:

tuple[float, float]

Notes

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

Source code at slepc4py/SLEPc/PEP.pyx:764

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:

int

Source code at slepc4py/SLEPc/PEP.pyx:1648

getJDFix()#

Get threshold for changing the target in the correction equation.

Not collective.

Returns:

The threshold for changing the target.

Return type:

float

See also

setJDFix, PEPJDGetFix

Source code at slepc4py/SLEPc/PEP.pyx:2566

getJDMinimalityIndex()#

Get the maximum allowed value of the minimality index.

Not collective.

Returns:

The maximum minimality index.

Return type:

int

Source code at slepc4py/SLEPc/PEP.pyx:2653

getJDProjection()#

Get the type of projection to be used in the Jacobi-Davidson solver.

Not collective.

Returns:

The type of projection.

Return type:

JDProjection

Source code at slepc4py/SLEPc/PEP.pyx:2690

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:

float

Source code at slepc4py/SLEPc/PEP.pyx:2523

getJDReusePreconditioner()#

Get the flag for reusing the preconditioner.

Not collective.

Returns:

The reuse flag.

Return type:

bool

Source code at slepc4py/SLEPc/PEP.pyx:2610

getLinearEPS()#

Get the eigensolver object associated to the polynomial eigenvalue solver.

Collective.

Returns:

The linear eigensolver.

Return type:

EPS

Source code at slepc4py/SLEPc/PEP.pyx:1912

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:

bool

Source code at slepc4py/SLEPc/PEP.pyx:1993

getLinearLinearization()#

Get the coeffs. defining the linearization of a quadratic eigenproblem.

Not collective.

Returns:

  • alpha (float) – First parameter of the linearization.

  • beta (float) – Second parameter of the linearization.

Return type:

tuple[float, float]

Source code at slepc4py/SLEPc/PEP.pyx:1953

getMonitor()#

Get the list of monitor functions.

Not collective.

Returns:

The list of monitor functions.

Return type:

PEPMonitorFunction

See also

setMonitor

Source code at slepc4py/SLEPc/PEP.pyx:1572

getOperators()#

Get the matrices associated with the eigenvalue problem.

Collective.

Returns:

The matrices associated with the eigensystem.

Return type:

list of petsc4py.PETSc.Mat

Source code at slepc4py/SLEPc/PEP.pyx:1361

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:

str

Source code at slepc4py/SLEPc/PEP.pyx:417

getProblemType()#

Get the problem type from the PEP object.

Not collective.

Returns:

The problem type that was previously set.

Return type:

ProblemType

Source code at slepc4py/SLEPc/PEP.pyx:553

getQArnoldiLocking()#

Get the locking flag used in the Q-Arnoldi method.

Not collective.

Returns:

The locking flag.

Return type:

bool

Source code at slepc4py/SLEPc/PEP.pyx:2084

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:

float

Source code at slepc4py/SLEPc/PEP.pyx:2040

getRG()#

Get the region object associated to the eigensolver.

Not collective.

Returns:

The region context.

Return type:

RG

See also

setRG, PEPGetRG

Source code at slepc4py/SLEPc/PEP.pyx:1287

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:

tuple[Refine, int, float, int, RefineScheme]

Source code at slepc4py/SLEPc/PEP.pyx:859

getRefineKSP()#

Get the KSP object used by the eigensolver in the refinement phase.

Collective.

Returns:

The linear solver object.

Return type:

petsc4py.PETSc.KSP

Source code at slepc4py/SLEPc/PEP.pyx:934

getST()#

Get the spectral transformation object associated to the eigensolver.

Not collective.

Returns:

The spectral transformation.

Return type:

ST

See also

setST, PEPGetST

Source code at slepc4py/SLEPc/PEP.pyx:1114

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:

bool

Source code at slepc4py/SLEPc/PEP.pyx:2476

getSTOARDetectZeros()#

Get the flag that enforces zero detection in spectrum slicing.

Not collective.

Returns:

The zero detection flag.

Return type:

bool

Source code at slepc4py/SLEPc/PEP.pyx:2316

getSTOARDimensions()#

Get the dimensions used for each subsolve step.

Not collective.

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.

Return type:

tuple[int, int, int]

Source code at slepc4py/SLEPc/PEP.pyx:2375

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:

  • shifts (ArrayReal) – The values of the shifts used internally in the solver.

  • inertias (ArrayInt) – The values of the inertia in each shift.

Return type:

tuple[ArrayReal, ArrayInt]

Notes

This call makes sense only for spectrum slicing runs, that is, when an interval has been given with setInterval() and SINVERT is set.

If called after solve(), all shifts used internally by the solver are returned (including both endpoints and any intermediate ones). If called before solve() and after setUp() then only the information of the endpoints of subintervals is available.

Source code at slepc4py/SLEPc/PEP.pyx:2400

getSTOARLinearization()#

Get the coefficients that define the linearization of a quadratic eigenproblem.

Not collective.

Returns:

  • alpha (float) – First parameter of the linearization.

  • beta (float) – Second parameter of the linearization.

Return type:

tuple[float, float]

Source code at slepc4py/SLEPc/PEP.pyx:2217

getSTOARLocking()#

Get the locking flag used in the STOAR method.

Not collective.

Returns:

The locking flag.

Return type:

bool

Source code at slepc4py/SLEPc/PEP.pyx:2264

getScale(Dl=None, Dr=None)#

Get the strategy used for scaling the polynomial eigenproblem.

Not collective.

Parameters:
Returns:

  • scale (Scale) – The scaling strategy.

  • alpha (float) – The scaling factor.

  • its (int) – The number of iterations of diagonal scaling.

  • lbda (float) – Approximation of the wanted eigenvalues (modulus).

Return type:

tuple[Scale, float, int, float]

See also

setScale, PEPGetScale

Source code at slepc4py/SLEPc/PEP.pyx:1151

getStoppingTest()#

Get the stopping test function.

Not collective.

Returns:

The stopping test function.

Return type:

PEPStoppingFunction

See also

setStoppingTest

Source code at slepc4py/SLEPc/PEP.pyx:1483

getTOARLocking()#

Get the locking flag used in the TOAR method.

Not collective.

Returns:

The locking flag.

Return type:

bool

Source code at slepc4py/SLEPc/PEP.pyx:2175

getTOARRestart()#

Get the restart parameter used in the TOAR method.

Not collective.

Returns:

The number of vectors to be kept at restart.

Return type:

float

Source code at slepc4py/SLEPc/PEP.pyx:2131

getTarget()#

Get the value of the target.

Not collective.

Returns:

The value of the target.

Return type:

Scalar

Notes

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

Source code at slepc4py/SLEPc/PEP.pyx:658

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:

  • tol (float) – The convergence tolerance.

  • max_it (int) – The maximum number of iterations.

Return type:

tuple[float, int]

Source code at slepc4py/SLEPc/PEP.pyx:708

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:

bool

Source code at slepc4py/SLEPc/PEP.pyx:998

getType()#

Get the PEP type of this object.

Not collective.

Returns:

The solver currently being used.

Return type:

str

See also

setType, PEPGetType

Source code at slepc4py/SLEPc/PEP.pyx:398

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:

Which

Source code at slepc4py/SLEPc/PEP.pyx:602

reset()#

Reset the PEP object.

Collective.

See also

PEPReset

Source code at slepc4py/SLEPc/PEP.pyx:338

Return type:

None

setBV(bv)#

Set a basis vectors object associated to the eigensolver.

Collective.

Parameters:

bv (BV) – The basis vectors context.

Return type:

None

See also

getBV, PEPSetBV

Source code at slepc4py/SLEPc/PEP.pyx:1270

setBasis(basis)#

Set the type of polynomial basis used.

Logically collective.

Set the type of polynomial basis used to describe the polynomial eigenvalue problem.

Parameters:

basis (Basis) – The basis to be set.

Return type:

None

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.

Source code at slepc4py/SLEPc/PEP.pyx:524

setCISSExtraction(extraction)#

Set the extraction technique used in the CISS solver.

Logically collective.

Parameters:

extraction (CISSExtraction) – The extraction technique.

Return type:

None

Source code at slepc4py/SLEPc/PEP.pyx:2711

setCISSRefinement(inner=None, blsize=None)#

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

Logically collective.

Parameters:
  • inner (int | None) – Number of iterative refinement iterations (inner loop).

  • blsize (int | None) – Number of iterative refinement iterations (blocksize loop).

Return type:

None

Source code at slepc4py/SLEPc/PEP.pyx:2881

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:
  • ip (int | None) – Number of integration points.

  • bs (int | None) – Block size.

  • ms (int | None) – Moment size.

  • npart (int | None) – Number of partitions when splitting the communicator.

  • bsmax (int | None) – Maximum block size.

  • realmats (bool) – True if A and B are real.

Return type:

None

Notes

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

Source code at slepc4py/SLEPc/PEP.pyx:2748

setCISSThreshold(delta=None, spur=None)#

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

Logically collective.

Parameters:
  • delta (float | None) – Threshold for numerical rank.

  • spur (float | None) – Spurious threshold (to discard spurious eigenpairs).

Return type:

None

Source code at slepc4py/SLEPc/PEP.pyx:2836

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:

None

Source code at slepc4py/SLEPc/PEP.pyx:840

setDS(ds)#

Set a direct solver object associated to the eigensolver.

Collective.

Parameters:

ds (DS) – The direct solver context.

Return type:

None

See also

getDS, PEPSetDS

Source code at slepc4py/SLEPc/PEP.pyx:1344

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

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

Logically collective.

Parameters:
  • nev (int | None) – Number of eigenvalues to compute.

  • ncv (int | None) – Maximum dimension of the subspace to be used by the solver.

  • mpd (int | None) – Maximum dimension allowed for the projected problem.

Return type:

None

Notes

Use DETERMINE 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.

Source code at slepc4py/SLEPc/PEP.pyx:1063

setEigenvalueComparison(comparison, args=None, kargs=None)#

Set an eigenvalue comparison function.

Logically collective.

Notes

This eigenvalue comparison function is used when setWhichEigenpairs() is set to PEP.Which.USER.

Source code at slepc4py/SLEPc/PEP.pyx:1500

Parameters:
Return type:

None

setExtract(extract)#

Set the extraction strategy to be used.

Logically collective.

Parameters:

extract (Extract) – The extraction strategy.

Return type:

None

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.

Source code at slepc4py/SLEPc/PEP.pyx:954

setFromOptions()#

Set PEP options from the options database.

Collective.

Notes

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

This routine must be called before setUp() if the user is to be allowed to set the solver type.

Source code at slepc4py/SLEPc/PEP.pyx:486

Return type:

None

setInitialSpace(space)#

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

Collective.

Parameters:

space (Vec | list[Vec]) – The initial space.

Return type:

None

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.

Source code at slepc4py/SLEPc/PEP.pyx:1420

setInterval(inta, intb)#

Set the computational interval for spectrum slicing.

Logically collective.

Parameters:
  • inta (float) – The left end of the interval.

  • intb (float) – The right end of the interval.

Return type:

None

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, see setWhichEigenpairs(). Note that in polynomial eigenproblems spectrum slicing is implemented in STOAR only.

Source code at slepc4py/SLEPc/PEP.pyx:790

setJDFix(fix)#

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

Logically collective.

Parameters:

fix (float) – Threshold for changing the target.

Return type:

None

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

getJDFix, PEPJDSetFix

Source code at slepc4py/SLEPc/PEP.pyx:2542

setJDMinimalityIndex(flag)#

Set the maximum allowed value for the minimality index.

Logically collective.

Parameters:

flag (int) – The maximum minimality index.

Return type:

None

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.

Source code at slepc4py/SLEPc/PEP.pyx:2629

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:

None

Source code at slepc4py/SLEPc/PEP.pyx:2672

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.

Parameters:

keep (float) – The number of vectors to be kept at restart.

Return type:

None

Notes

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

Source code at slepc4py/SLEPc/PEP.pyx:2497

setJDReusePreconditioner(flag)#

Set a flag indicating whether the preconditioner must be reused or not.

Logically collective.

Parameters:

flag (bool) – The reuse flag.

Return type:

None

Notes

The default value is False. If set to True, the preconditioner is built only at the beginning, using the target value. Otherwise, it may be rebuilt (depending on the fix parameter) at each iteration from the Ritz value.

Source code at slepc4py/SLEPc/PEP.pyx:2585

setLinearEPS(eps)#

Set an eigensolver object associated to the polynomial eigenvalue solver.

Collective.

Parameters:

eps (EPS) – The linear eigensolver.

Return type:

None

Source code at slepc4py/SLEPc/PEP.pyx:1895

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:

None

Source code at slepc4py/SLEPc/PEP.pyx:1975

setLinearLinearization(alpha=1.0, beta=0.0)#

Set the coefficients that define the linearization of a quadratic eigenproblem.

Logically collective.

Parameters:
  • alpha (float) – First parameter of the linearization.

  • beta (float) – Second parameter of the linearization.

Return type:

None

Source code at slepc4py/SLEPc/PEP.pyx:1932

setMonitor(monitor, args=None, kargs=None)#

Append a monitor function to the list of monitors.

Logically collective.

Source code at slepc4py/SLEPc/PEP.pyx:1547

Parameters:
Return type:

None

setOperators(operators)#

Set the matrices associated with the eigenvalue problem.

Collective.

Parameters:

operators (list[Mat]) – The matrices associated with the eigensystem.

Return type:

None

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.

Source code at slepc4py/SLEPc/PEP.pyx:1387

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:

None

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_")

Source code at slepc4py/SLEPc/PEP.pyx:436

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:

None

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.

Source code at slepc4py/SLEPc/PEP.pyx:572

setQArnoldiLocking(lock)#

Toggle between locking and non-locking variants of the Q-Arnoldi method.

Logically collective.

Parameters:

lock (bool) – True if the locking variant must be selected.

Return type:

None

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

Source code at slepc4py/SLEPc/PEP.pyx:2059

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.

Parameters:

keep (float) – The number of vectors to be kept at restart.

Return type:

None

Notes

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

Source code at slepc4py/SLEPc/PEP.pyx:2014

setRG(rg)#

Set a region object associated to the eigensolver.

Collective.

Parameters:

rg (RG) – The region context.

Return type:

None

See also

getRG, PEPSetRG

Source code at slepc4py/SLEPc/PEP.pyx:1307

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:
  • ref (Refine) – The refinement type.

  • npart (int | None) – The number of partitions of the communicator.

  • tol (float | None) – The convergence tolerance.

  • its (int | None) – The maximum number of refinement iterations.

  • scheme (RefineScheme | None) – Scheme for solving linear systems.

Return type:

None

Source code at slepc4py/SLEPc/PEP.pyx:890

setST(st)#

Set a spectral transformation object associated to the eigensolver.

Collective.

Parameters:

st (ST) – The spectral transformation.

Return type:

None

See also

getST, PEPSetST

Source code at slepc4py/SLEPc/PEP.pyx:1134

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.

Parameters:

flag (bool) – Whether the eigenvalue type is checked or not.

Return type:

None

Notes

This option is relevant only for spectrum slicing computations, but is ignored in slepc4py.SLEPc.PEP.ProblemType.HYPERBOLIC problems.

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.

Source code at slepc4py/SLEPc/PEP.pyx:2444

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.

Parameters:

detect (bool) – True if zeros must checked for.

Return type:

None

Notes

This call makes sense only for spectrum slicing runs, that is, when an interval has been given with setInterval() and SINVERT is 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.

Source code at slepc4py/SLEPc/PEP.pyx:2283

setSTOARDimensions(nev=None, ncv=None, mpd=None)#

Set the dimensions used for each subsolve step.

Logically collective.

Parameters:
  • nev (int | None) – Number of eigenvalues to compute.

  • ncv (int | None) – Maximum dimension of the subspace to be used by the solver.

  • mpd (int | None) – Maximum dimension allowed for the projected problem.

Return type:

None

Notes

This call makes sense only for spectrum slicing runs, that is, when an interval has been given with setInterval() and SINVERT is set.

The meaning of the parameters is the same as in setDimensions(), but the ones here apply to every subsolve done by the child PEP object.

Source code at slepc4py/SLEPc/PEP.pyx:2335

setSTOARLinearization(alpha=1.0, beta=0.0)#

Set the coefficients that define the linearization of a quadratic eigenproblem.

Logically collective.

Parameters:
  • alpha (float) – First parameter of the linearization.

  • beta (float) – Second parameter of the linearization.

Return type:

None

Source code at slepc4py/SLEPc/PEP.pyx:2196

setSTOARLocking(lock)#

Toggle between locking and non-locking variants of the STOAR method.

Logically collective.

Parameters:

lock (bool) – True if the locking variant must be selected.

Return type:

None

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

Source code at slepc4py/SLEPc/PEP.pyx:2239

setScale(scale, alpha=None, Dl=None, Dr=None, its=None, lbda=None)#

Set the scaling strategy to be used.

Collective.

Parameters:
Return type:

None

See also

getScale, PEPSetScale

Source code at slepc4py/SLEPc/PEP.pyx:1204

setStoppingTest(stopping, args=None, kargs=None)#

Set a function to decide when to stop the outer iteration of the eigensolver.

Logically collective.

Source code at slepc4py/SLEPc/PEP.pyx:1459

Parameters:
Return type:

None

setTOARLocking(lock)#

Toggle between locking and non-locking variants of the TOAR method.

Logically collective.

Parameters:

lock (bool) – True if the locking variant must be selected.

Return type:

None

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

Source code at slepc4py/SLEPc/PEP.pyx:2150

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.

Parameters:

keep (float) – The number of vectors to be kept at restart.

Return type:

None

Notes

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

Source code at slepc4py/SLEPc/PEP.pyx:2105

setTarget(target)#

Set the value of the target.

Logically collective.

Parameters:

target (Scalar) – The value of the target.

Return type:

None

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.

Source code at slepc4py/SLEPc/PEP.pyx:681

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:
  • tol (float | None) – The convergence tolerance.

  • max_it (int | None) – The maximum number of iterations

Return type:

None

Notes

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

Source code at slepc4py/SLEPc/PEP.pyx:733

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.

Parameters:

trackall (bool) – Whether to compute all residuals or not.

Return type:

None

Source code at slepc4py/SLEPc/PEP.pyx:1017

setType(pep_type)#

Set the particular solver to be used in the PEP object.

Logically collective.

Parameters:

pep_type (Type | str) – The solver to be used.

Return type:

None

Notes

The default is TOAR. Normally, it is best to use setFromOptions() 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

getType, PEPSetType

Source code at slepc4py/SLEPc/PEP.pyx:371

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 ST object.

This function need not be called explicitly in most cases, since solve() calls it. It can be useful when one wants to measure the set-up time separately from the solve time.

See also

solve, PEPSetUp

Source code at slepc4py/SLEPc/PEP.pyx:1604

Return type:

None

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:

None

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

The target is a scalar value provided with setTarget().

The criterion PEP.Which.TARGET_IMAGINARY is available only in case PETSc and SLEPc have been built with complex scalars.

PEP.Which.ALL is intended for use in combination with an interval (see setInterval()), when all eigenvalues within the interval are requested, or in the context of the PEP.Type.CISS solver for computing all eigenvalues in a region.

Source code at slepc4py/SLEPc/PEP.pyx:621

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. Call getConverged() to get the actual number of computed solutions, and getConvergedReason() to determine if the solver converged or failed and why.

Source code at slepc4py/SLEPc/PEP.pyx:1626

Return type:

None

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:

None

Source code at slepc4py/SLEPc/PEP.pyx:1855

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:

None

Source code at slepc4py/SLEPc/PEP.pyx:1874

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:

None

See also

PEPView

Source code at slepc4py/SLEPc/PEP.pyx:305

Attributes Documentation

bv#

The basis vectors (BV) object associated.

Source code at slepc4py/SLEPc/PEP.pyx:3011

ds#

The direct solver (DS) object associated.

Source code at slepc4py/SLEPc/PEP.pyx:3025

extract#

The type of extraction technique to be employed.

Source code at slepc4py/SLEPc/PEP.pyx:2976

max_it#

The maximum iteration count.

Source code at slepc4py/SLEPc/PEP.pyx:2990

problem_type#

The type of the eigenvalue problem.

Source code at slepc4py/SLEPc/PEP.pyx:2955

rg#

The region (RG) object associated.

Source code at slepc4py/SLEPc/PEP.pyx:3018

st#

The spectral transformation (ST) object associated.

Source code at slepc4py/SLEPc/PEP.pyx:3004

target#

The value of the target.

Source code at slepc4py/SLEPc/PEP.pyx:2969

tol#

The tolerance.

Source code at slepc4py/SLEPc/PEP.pyx:2983

track_all#

Compute the residual norm of all approximate eigenpairs.

Source code at slepc4py/SLEPc/PEP.pyx:2997

which#

The portion of the spectrum to be sought.

Source code at slepc4py/SLEPc/PEP.pyx:2962