slepc4py.SLEPc.NEP#

class slepc4py.SLEPc.NEP#

Bases: Object

Nonlinear Eigenvalue Problem Solver.

The Nonlinear Eigenvalue Problem (NEP) solver is the object provided by slepc4py for specifying an eigenvalue problem that is nonlinear with respect to the eigenvalue (not the eigenvector). This is intended for general nonlinear problems (rather than polynomial eigenproblems) described as \(T(\lambda) x=0\).

Enumerations

CISSExtraction

NEP CISS extraction technique.

Conv

NEP convergence test.

ConvergedReason

NEP convergence reasons.

ErrorType

NEP error type to assess accuracy of computed solutions.

ProblemType

NEP problem type.

Refine

NEP refinement strategy.

RefineScheme

NEP scheme for solving linear systems during iterative refinement.

Stop

NEP stopping test.

Type

NEP type.

Which

NEP desired part of spectrum.

Methods Summary

appendOptionsPrefix([prefix])

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

applyResolvent(omega, v, r[, rg])

Apply the resolvent \(T^{-1}(z)\) to a given vector.

cancelMonitor()

Clear all monitors for a NEP object.

computeError(i[, etype])

Compute the error associated with the i-th computed eigenpair.

create([comm])

Create the NEP object.

destroy()

Destroy the NEP object.

errorView([etype, viewer])

Display the errors associated with the computed solution.

getBV()

Get the basis vectors object associated to the eigensolver.

getCISSExtraction()

Get the extraction technique used in the CISS solver.

getCISSKSPs()

Get the list 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.

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.

getFunction()

Get the function to compute the nonlinear Function \(T(\lambda)\).

getInterpolInterpolation()

Get the tolerance and maximum degree for the interpolation polynomial.

getInterpolPEP()

Get the associated polynomial eigensolver object.

getIterationNumber()

Get the current iteration number.

getJacobian()

Get the function to compute the Jacobian \(T'(\lambda)\) and J.

getLeftEigenvector(i, Wr[, Wi])

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

getMonitor()

Get the list of monitor functions.

getNArnoldiKSP()

Get the linear solver object associated with the nonlinear eigensolver.

getNArnoldiLagPreconditioner()

Get how often the preconditioner is rebuilt.

getNLEIGSEPS()

Get the linear eigensolver object associated with the nonlinear eigensolver.

getNLEIGSFullBasis()

Get the flag that indicates if NLEIGS is using the full-basis variant.

getNLEIGSInterpolation()

Get the tolerance and maximum degree for the interpolation polynomial.

getNLEIGSKSPs()

Get the list of linear solver objects associated with the NLEIGS solver.

getNLEIGSLocking()

Get the locking flag used in the NLEIGS method.

getNLEIGSRKShifts()

Get the list of shifts used in the Rational Krylov method.

getNLEIGSRestart()

Get the restart parameter used in the NLEIGS method.

getOptionsPrefix()

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

getProblemType()

Get the problem type from the NEP object.

getRG()

Get the region object associated to the eigensolver.

getRIIConstCorrectionTol()

Get the constant tolerance flag.

getRIIDeflationThreshold()

Get the threshold value that controls deflation.

getRIIHermitian()

Get if the Hermitian version must be used by the solver.

getRIIKSP()

Get the linear solver object associated with the nonlinear eigensolver.

getRIILagPreconditioner()

Get how often the preconditioner is rebuilt.

getRIIMaximumIterations()

Get the maximum number of inner iterations of RII.

getRefine()

Get the refinement strategy used by the NEP object.

getRefineKSP()

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

getSLPDeflationThreshold()

Get the threshold value that controls deflation.

getSLPEPS()

Get the linear eigensolver object associated with the nonlinear eigensolver.

getSLPEPSLeft()

Get the left eigensolver.

getSLPKSP()

Get the linear solver object associated with the nonlinear eigensolver.

getSplitOperator()

Get the operator of the nonlinear eigenvalue problem in split form.

getSplitPreconditioner()

Get the operator of the split preconditioner.

getStoppingTest()

Get the stopping test function.

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.

getTwoSided()

Get the flag indicating if a two-sided variant is being used.

getType()

Get the NEP type of this object.

getWhichEigenpairs()

Get which portion of the spectrum is to be sought.

reset()

Reset the NEP object.

setBV(bv)

Set the basis vectors object associated to the eigensolver.

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.

setEigenvalueComparison(comparison[, args, ...])

Set an eigenvalue comparison function.

setFromOptions()

Set NEP options from the options database.

setFunction(function[, F, P, args, kargs])

Set the function to compute the nonlinear Function \(T(\lambda)\).

setInitialSpace(space)

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

setInterpolInterpolation([tol, deg])

Set the tolerance and maximum degree for the interpolation polynomial.

setInterpolPEP(pep)

Set a polynomial eigensolver object associated to the nonlinear eigensolver.

setJacobian(jacobian[, J, args, kargs])

Set the function to compute the Jacobian \(T'(\lambda)\).

setMonitor(monitor[, args, kargs])

Append a monitor function to the list of monitors.

setNArnoldiKSP(ksp)

Set a linear solver object associated to the nonlinear eigensolver.

setNArnoldiLagPreconditioner(lag)

Set when the preconditioner is rebuilt in the nonlinear solve.

setNLEIGSEPS(eps)

Set a linear eigensolver object associated to the nonlinear eigensolver.

setNLEIGSFullBasis([fullbasis])

Set TOAR-basis (default) or full-basis variants of the NLEIGS method.

setNLEIGSInterpolation([tol, deg])

Set the tolerance and maximum degree for the interpolation polynomial.

setNLEIGSLocking(lock)

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

setNLEIGSRKShifts(shifts)

Set a list of shifts to be used in the Rational Krylov method.

setNLEIGSRestart(keep)

Set the restart parameter for the NLEIGS method.

setOptionsPrefix([prefix])

Set the prefix used for searching for all NEP options in the database.

setProblemType(problem_type)

Set the type of the eigenvalue problem.

setRG(rg)

Set a region object associated to the eigensolver.

setRIIConstCorrectionTol(cct)

Set a flag to keep the tolerance used in the linear solver constant.

setRIIDeflationThreshold(deftol)

Set the threshold used to switch between deflated and non-deflated.

setRIIHermitian(herm)

Set a flag to use the Hermitian version of the solver.

setRIIKSP(ksp)

Set a linear solver object associated to the nonlinear eigensolver.

setRIILagPreconditioner(lag)

Set when the preconditioner is rebuilt in the nonlinear solve.

setRIIMaximumIterations(its)

Set the max.

setRefine(ref[, npart, tol, its, scheme])

Set the refinement strategy used by the NEP object.

setSLPDeflationThreshold(deftol)

Set the threshold used to switch between deflated and non-deflated.

setSLPEPS(eps)

Set a linear eigensolver object associated to the nonlinear eigensolver.

setSLPEPSLeft(eps)

Set a linear eigensolver object associated to the nonlinear eigensolver.

setSLPKSP(ksp)

Set a linear solver object associated to the nonlinear eigensolver.

setSplitOperator(A, f[, structure])

Set the operator of the nonlinear eigenvalue problem in split form.

setSplitPreconditioner(P[, structure])

Set the operator in split form.

setStoppingTest(stopping[, args, kargs])

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

setTarget(target)

Set the value of the target.

setTolerances([tol, maxit])

Set the tolerance and max.

setTrackAll(trackall)

Set if the solver must compute the residual of all approximate eigenpairs.

setTwoSided(twosided)

Set the solver to use a two-sided variant.

setType(nep_type)

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

setUp()

Set up all the internal data structures.

setWhichEigenpairs(which)

Set which portion of the spectrum is to be sought.

solve()

Solve the nonlinear eigenproblem.

valuesView([viewer])

Display the computed eigenvalues in a viewer.

vectorsView([viewer])

Output computed eigenvectors to a viewer.

view([viewer])

Print the NEP data structure.

Attributes Summary

bv

The basis vectors (BV) object associated.

ds

The direct solver (DS) object associated.

max_it

The maximum iteration count used by the NEP convergence tests.

problem_type

The problem type from the NEP object.

rg

The region (RG) object associated.

target

The value of the target.

tol

The tolerance used by the NEP convergence tests.

track_all

Compute the residual 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 NEP options in the database.

Logically collective.

Parameters:

prefix (str | None) – The prefix string to prepend to all NEP option requests.

Return type:

None

Source code at slepc4py/SLEPc/NEP.pyx:385

applyResolvent(omega, v, r, rg=None)#

Apply the resolvent \(T^{-1}(z)\) to a given vector.

Collective.

Parameters:
  • omega (Scalar) – Value where the resolvent must be evaluated.

  • v (Vec) – Input vector.

  • r (Vec) – Placeholder for the result vector.

  • rg (RG | None) – Region.

Return type:

None

Notes

The resolvent \(T^{-1}(z)=\sum_i(z-\lambda_i)^{-1}x_iy_i^*\) is evaluated at \(z=\omega\) and the matrix-vector product \(r = T^{-1}(\omega) v\) is computed. Vectors \(x_i,y_i\) are right and left eigenvectors, respectively, normalized so that \(y_i^*T'(\lambda_i)x_i=1\). The sum contains only eigenvectors that have been previously computed with solve(), and if a region rg is given then only those corresponding to eigenvalues inside the region are considered.

Source code at slepc4py/SLEPc/NEP.pyx:1834

cancelMonitor()#

Clear all monitors for a NEP object.

Logically collective.

See also

NEPMonitorCancel

Source code at slepc4py/SLEPc/NEP.pyx:1159

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

If the computation of left eigenvectors was enabled with setTwoSided(), then the error will be computed using the maximum of the value above and the left residual norm \(\|y^*T(\lambda)\|_2\), where \(y\) is the approximate left eigenvector.

Source code at slepc4py/SLEPc/NEP.pyx:1385

create(comm=None)#

Create the NEP object.

Collective.

Parameters:

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

Return type:

Self

See also

NEPCreate

Source code at slepc4py/SLEPc/NEP.pyx:268

destroy()#

Destroy the NEP object.

Collective.

See also

NEPDestroy

Source code at slepc4py/SLEPc/NEP.pyx:242

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/NEP.pyx:1428

getBV()#

Get the basis vectors object associated to the eigensolver.

Not collective.

Returns:

The basis vectors context.

Return type:

BV

See also

setBV, NEPGetBV

Source code at slepc4py/SLEPc/NEP.pyx:877

getCISSExtraction()#

Get the extraction technique used in the CISS solver.

Not collective.

Returns:

The extraction technique.

Return type:

CISSExtraction

Source code at slepc4py/SLEPc/NEP.pyx:2792

getCISSKSPs()#

Get the list 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/NEP.pyx:2989

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/NEP.pyx:2967

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/NEP.pyx:2865

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/NEP.pyx:2922

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/NEP.pyx:1256

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/NEP.pyx:1237

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/NEP.pyx:624

getDS()#

Get the direct solver associated to the eigensolver.

Not collective.

Returns:

The direct solver context.

Return type:

DS

See also

setDS, NEPGetDS

Source code at slepc4py/SLEPc/NEP.pyx:951

getDimensions()#

Get the number of eigenvalues to compute.

Not collective.

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

Return type:

tuple[int, int, int]

Source code at slepc4py/SLEPc/NEP.pyx:795

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/NEP.pyx:1282

getEigenvalueComparison()#

Get the eigenvalue comparison function.

Not collective.

Returns:

The eigenvalue comparison function.

Return type:

NEPEigenvalueComparison

Source code at slepc4py/SLEPc/NEP.pyx:1100

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/NEP.pyx:1356

getFunction()#

Get the function to compute the nonlinear Function \(T(\lambda)\).

Collective.

Get the function to compute the nonlinear Function \(T(\lambda)\) and the matrix.

Returns:

Return type:

tuple[petsc4py.PETSc.Mat, petsc4py.PETSc.Mat, NEPFunction]

Source code at slepc4py/SLEPc/NEP.pyx:1540

getInterpolInterpolation()#

Get the tolerance and maximum degree for the interpolation polynomial.

Not collective.

Returns:

  • tol (float) – The tolerance to stop computing polynomial coefficients.

  • deg (int) – The maximum degree of interpolation.

Return type:

tuple[float, int]

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

getInterpolPEP()#

Get the associated polynomial eigensolver object.

Collective.

Returns:

The polynomial eigensolver.

Return type:

PEP

Source code at slepc4py/SLEPc/NEP.pyx:2401

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/NEP.pyx:1215

getJacobian()#

Get the function to compute the Jacobian \(T'(\lambda)\) and J.

Collective.

Get the function to compute the Jacobian \(T'(\lambda)\) and the matrix.

Returns:

Return type:

tuple[petsc4py.PETSc.Mat, NEPJacobian]

Source code at slepc4py/SLEPc/NEP.pyx:1606

getLeftEigenvector(i, Wr, Wi=None)#

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

Collective.

Parameters:
  • i (int) – Index of the solution to be obtained.

  • Wr (Vec) – Placeholder for the returned eigenvector (real part).

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

Return type:

None

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

Source code at slepc4py/SLEPc/NEP.pyx:1323

getMonitor()#

Get the list of monitor functions.

Not collective.

Returns:

The list of monitor functions.

Return type:

NEPMonitorFunction

See also

setMonitor

Source code at slepc4py/SLEPc/NEP.pyx:1142

getNArnoldiKSP()#

Get the linear solver object associated with the nonlinear eigensolver.

Collective.

Returns:

The linear solver object.

Return type:

petsc4py.PETSc.KSP

Source code at slepc4py/SLEPc/NEP.pyx:2318

getNArnoldiLagPreconditioner()#

Get how often the preconditioner is rebuilt.

Not collective.

Returns:

The lag parameter.

Return type:

int

Source code at slepc4py/SLEPc/NEP.pyx:2363

getNLEIGSEPS()#

Get the linear eigensolver object associated with the nonlinear eigensolver.

Collective.

Returns:

The linear eigensolver.

Return type:

EPS

Source code at slepc4py/SLEPc/NEP.pyx:2675

getNLEIGSFullBasis()#

Get the flag that indicates if NLEIGS is using the full-basis variant.

Not collective.

Returns:

True if the full-basis variant is selected.

Return type:

bool

Source code at slepc4py/SLEPc/NEP.pyx:2639

getNLEIGSInterpolation()#

Get the tolerance and maximum degree for the interpolation polynomial.

Not collective.

Get the tolerance and maximum degree when building the interpolation via divided differences.

Returns:

  • tol (float) – The tolerance to stop computing divided differences.

  • deg (int) – The maximum degree of interpolation.

Return type:

tuple[float, int]

Source code at slepc4py/SLEPc/NEP.pyx:2581

getNLEIGSKSPs()#

Get the list of linear solver objects associated with the NLEIGS 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 shifts provided by the user, or 1 if the user did not provide shifts.

Source code at slepc4py/SLEPc/NEP.pyx:2747

getNLEIGSLocking()#

Get the locking flag used in the NLEIGS method.

Not collective.

Returns:

The locking flag.

Return type:

bool

Source code at slepc4py/SLEPc/NEP.pyx:2536

getNLEIGSRKShifts()#

Get the list of shifts used in the Rational Krylov method.

Not collective.

Returns:

The shift values.

Return type:

ArrayScalar

Source code at slepc4py/SLEPc/NEP.pyx:2722

getNLEIGSRestart()#

Get the restart parameter used in the NLEIGS method.

Not collective.

Returns:

The number of vectors to be kept at restart.

Return type:

float

Source code at slepc4py/SLEPc/NEP.pyx:2492

getOptionsPrefix()#

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

Not collective.

Returns:

The prefix string set for this NEP object.

Return type:

str

Source code at slepc4py/SLEPc/NEP.pyx:335

getProblemType()#

Get the problem type from the NEP object.

Not collective.

Returns:

The problem type that was previously set.

Return type:

ProblemType

Source code at slepc4py/SLEPc/NEP.pyx:423

getRG()#

Get the region object associated to the eigensolver.

Not collective.

Returns:

The region context.

Return type:

RG

See also

setRG, NEPGetRG

Source code at slepc4py/SLEPc/NEP.pyx:914

getRIIConstCorrectionTol()#

Get the constant tolerance flag.

Not collective.

Returns:

If True, the petsc4py.PETSc.KSP relative tolerance is constant.

Return type:

bool

Source code at slepc4py/SLEPc/NEP.pyx:1942

getRIIDeflationThreshold()#

Get the threshold value that controls deflation.

Not collective.

Returns:

The threshold value.

Return type:

float

Source code at slepc4py/SLEPc/NEP.pyx:2081

getRIIHermitian()#

Get if the Hermitian version must be used by the solver.

Not collective.

Returns:

If True, the Hermitian version is used.

Return type:

bool

Source code at slepc4py/SLEPc/NEP.pyx:2031

getRIIKSP()#

Get the linear solver object associated with the nonlinear eigensolver.

Collective.

Returns:

The linear solver object.

Return type:

petsc4py.PETSc.KSP

Source code at slepc4py/SLEPc/NEP.pyx:2117

getRIILagPreconditioner()#

Get how often the preconditioner is rebuilt.

Not collective.

Returns:

The lag parameter.

Return type:

int

Source code at slepc4py/SLEPc/NEP.pyx:1898

getRIIMaximumIterations()#

Get the maximum number of inner iterations of RII.

Not collective.

Returns:

Maximum inner iterations.

Return type:

int

Source code at slepc4py/SLEPc/NEP.pyx:1983

getRefine()#

Get the refinement strategy used by the NEP 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/NEP.pyx:663

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/NEP.pyx:738

getSLPDeflationThreshold()#

Get the threshold value that controls deflation.

Not collective.

Returns:

The threshold value.

Return type:

float

Source code at slepc4py/SLEPc/NEP.pyx:2167

getSLPEPS()#

Get the linear eigensolver object associated with the nonlinear eigensolver.

Collective.

Returns:

The linear eigensolver.

Return type:

EPS

Source code at slepc4py/SLEPc/NEP.pyx:2203

getSLPEPSLeft()#

Get the left eigensolver.

Collective.

Returns:

The linear eigensolver.

Return type:

EPS

Source code at slepc4py/SLEPc/NEP.pyx:2242

getSLPKSP()#

Get the linear solver object associated with the nonlinear eigensolver.

Collective.

Returns:

The linear solver object.

Return type:

petsc4py.PETSc.KSP

Source code at slepc4py/SLEPc/NEP.pyx:2279

getSplitOperator()#

Get the operator of the nonlinear eigenvalue problem in split form.

Collective.

Returns:

Return type:

tuple[list[petsc4py.PETSc.Mat], list[FN], petsc4py.PETSc.Mat.Structure]

Source code at slepc4py/SLEPc/NEP.pyx:1683

getSplitPreconditioner()#

Get the operator of the split preconditioner.

Not collective.

Returns:

Return type:

tuple[list[petsc4py.PETSc.Mat], petsc4py.PETSc.Mat.Structure]

Source code at slepc4py/SLEPc/NEP.pyx:1752

getStoppingTest()#

Get the stopping test function.

Not collective.

Returns:

The stopping test function.

Return type:

NEPStoppingFunction

See also

setStoppingTest

Source code at slepc4py/SLEPc/NEP.pyx:1053

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/NEP.pyx:521

getTolerances()#

Get the tolerance and maximum iteration count.

Not collective.

Get the tolerance and maximum iteration count used by the default NEP convergence tests.

Returns:

  • tol (float) – The convergence tolerance.

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

Return type:

tuple[float, int]

Source code at slepc4py/SLEPc/NEP.pyx:571

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/NEP.pyx:758

getTwoSided()#

Get the flag indicating if a two-sided variant is being used.

Not collective.

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

Returns:

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

Return type:

bool

Source code at slepc4py/SLEPc/NEP.pyx:1781

getType()#

Get the NEP type of this object.

Not collective.

Returns:

The solver currently being used.

Return type:

str

See also

setType, NEPGetType

Source code at slepc4py/SLEPc/NEP.pyx:316

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/NEP.pyx:467

reset()#

Reset the NEP object.

Collective.

See also

NEPReset

Source code at slepc4py/SLEPc/NEP.pyx:256

Return type:

None

setBV(bv)#

Set the basis vectors object associated to the eigensolver.

Collective.

Parameters:

bv (BV) – The basis vectors context.

Return type:

None

See also

getBV, NEPSetBV

Source code at slepc4py/SLEPc/NEP.pyx:897

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/NEP.pyx:2774

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/NEP.pyx:2944

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 NEP communicator. Otherwise, the communicator is split into npart communicators, so that npart petsc4py.PETSc.KSP solves proceed simultaneously.

Source code at slepc4py/SLEPc/NEP.pyx:2811

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/NEP.pyx:2899

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/NEP.pyx:644

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, NEPSetDS

Source code at slepc4py/SLEPc/NEP.pyx:971

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

Set the number of eigenvalues to compute.

Logically collective.

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

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/NEP.pyx:823

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 NEP.Which.USER.

Source code at slepc4py/SLEPc/NEP.pyx:1070

Parameters:
Return type:

None

setFromOptions()#

Set NEP 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/NEP.pyx:404

Return type:

None

setFunction(function, F=None, P=None, args=None, kargs=None)#

Set the function to compute the nonlinear Function \(T(\lambda)\).

Collective.

Set the function to compute the nonlinear Function \(T(\lambda)\) as well as the location to store the matrix.

Parameters:
Return type:

None

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

setInitialSpace(space)#

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

Collective.

Parameters:

space (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/NEP.pyx:990

setInterpolInterpolation(tol=None, deg=None)#

Set the tolerance and maximum degree for the interpolation polynomial.

Collective.

Parameters:
  • tol (float | None) – The tolerance to stop computing polynomial coefficients.

  • deg (int | None) – The maximum degree of interpolation.

Return type:

None

Source code at slepc4py/SLEPc/NEP.pyx:2421

setInterpolPEP(pep)#

Set a polynomial eigensolver object associated to the nonlinear eigensolver.

Collective.

Parameters:

pep (PEP) – The polynomial eigensolver.

Return type:

None

Source code at slepc4py/SLEPc/NEP.pyx:2384

setJacobian(jacobian, J=None, args=None, kargs=None)#

Set the function to compute the Jacobian \(T'(\lambda)\).

Collective.

Set the function to compute the Jacobian \(T'(\lambda)\) as well as the location to store the matrix.

Parameters:
Return type:

None

Source code at slepc4py/SLEPc/NEP.pyx:1570

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

Append a monitor function to the list of monitors.

Logically collective.

Source code at slepc4py/SLEPc/NEP.pyx:1117

Parameters:
Return type:

None

setNArnoldiKSP(ksp)#

Set a linear solver object associated to the nonlinear eigensolver.

Collective.

Parameters:

ksp (KSP) – The linear solver object.

Return type:

None

Source code at slepc4py/SLEPc/NEP.pyx:2301

setNArnoldiLagPreconditioner(lag)#

Set when the preconditioner is rebuilt in the nonlinear solve.

Logically collective.

Parameters:

lag (int) – 0 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within the nonlinear iteration, 2 means every second time the Jacobian is built, etc.

Return type:

None

Notes

The default is 1. The preconditioner is ALWAYS built in the first iteration of a nonlinear solve.

Source code at slepc4py/SLEPc/NEP.pyx:2338

setNLEIGSEPS(eps)#

Set a linear eigensolver object associated to the nonlinear eigensolver.

Collective.

Parameters:

eps (EPS) – The linear eigensolver.

Return type:

None

Source code at slepc4py/SLEPc/NEP.pyx:2658

setNLEIGSFullBasis(fullbasis=True)#

Set TOAR-basis (default) or full-basis variants of the NLEIGS method.

Logically collective.

Toggle between TOAR-basis (default) and full-basis variants of the NLEIGS method.

Parameters:

fullbasis (bool) – True if the full-basis variant must be selected.

Return type:

None

Notes

The default is to use a compact representation of the Krylov basis, that is, \(V = (I \otimes U) S\), with a BV of type TENSOR. This behavior can be changed so that the full basis \(V\) is explicitly stored and operated with. This variant is more expensive in terms of memory and computation, but is necessary in some cases, particularly for two-sided computations, see setTwoSided().

In the full-basis variant, the NLEIGS solver uses an EPS object to explicitly solve the linearized eigenproblem, see getNLEIGSEPS().

Source code at slepc4py/SLEPc/NEP.pyx:2606

setNLEIGSInterpolation(tol=None, deg=None)#

Set the tolerance and maximum degree for the interpolation polynomial.

Collective.

Set the tolerance and maximum degree when building the interpolation via divided differences.

Parameters:
  • tol (float | None) – The tolerance to stop computing divided differences.

  • deg (int | None) – The maximum degree of interpolation.

Return type:

None

Source code at slepc4py/SLEPc/NEP.pyx:2555

setNLEIGSLocking(lock)#

Toggle between locking and non-locking variants of the NLEIGS 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/NEP.pyx:2511

setNLEIGSRKShifts(shifts)#

Set a list of shifts to be used in the Rational Krylov method.

Collective.

Parameters:

shifts (Sequence[Scalar]) – Values specifying the shifts.

Return type:

None

Notes

If only one shift is provided, the built subspace is equivalent to shift-and-invert Krylov-Schur (provided that the absolute convergence criterion is used). Otherwise, the rational Krylov variant is run.

Source code at slepc4py/SLEPc/NEP.pyx:2695

setNLEIGSRestart(keep)#

Set the restart parameter for the NLEIGS method.

Logically collective.

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/NEP.pyx:2468

setOptionsPrefix(prefix=None)#

Set the prefix used for searching for all NEP options in the database.

Logically collective.

Parameters:

prefix (str | None) – The prefix string to prepend to all NEP 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 NEP contexts, one could call:

N1.setOptionsPrefix("nep1_")
N2.setOptionsPrefix("nep2_")

Source code at slepc4py/SLEPc/NEP.pyx:354

setProblemType(problem_type)#

Set the type of the eigenvalue problem.

Logically collective.

Parameters:

problem_type (ProblemType) – The problem type to be set.

Return type:

None

Notes

This function is used to provide a hint to the NEP solver to exploit certain properties of the nonlinear eigenproblem. This hint may be used or not, depending on the solver. By default, no particular structure is assumed.

Source code at slepc4py/SLEPc/NEP.pyx:442

setRG(rg)#

Set a region object associated to the eigensolver.

Collective.

Parameters:

rg (RG) – The region context.

Return type:

None

See also

getRG, NEPSetRG

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

setRIIConstCorrectionTol(cct)#

Set a flag to keep the tolerance used in the linear solver constant.

Logically collective.

Parameters:

cct (bool) – If True, the petsc4py.PETSc.KSP relative tolerance is constant.

Return type:

None

Notes

By default, an exponentially decreasing tolerance is set in the KSP used within the nonlinear iteration, so that each Newton iteration requests better accuracy than the previous one. The constant correction tolerance flag stops this behavior.

Source code at slepc4py/SLEPc/NEP.pyx:1917

setRIIDeflationThreshold(deftol)#

Set the threshold used to switch between deflated and non-deflated.

Logically collective.

Set the threshold value used to switch between deflated and non-deflated iteration.

Parameters:

deftol (float) – The threshold value.

Return type:

None

Notes

Normally, the solver iterates on the extended problem in order to deflate previously converged eigenpairs. If this threshold is set to a nonzero value, then once the residual error is below this threshold the solver will continue the iteration without deflation. The intention is to be able to improve the current eigenpair further, despite having previous eigenpairs with somewhat bad precision.

Source code at slepc4py/SLEPc/NEP.pyx:2050

setRIIHermitian(herm)#

Set a flag to use the Hermitian version of the solver.

Logically collective.

Set a flag to indicate if the Hermitian version of the scalar nonlinear equation must be used by the solver.

Parameters:

herm (bool) – If True, the Hermitian version is used.

Return type:

None

Notes

By default, the scalar nonlinear equation \(x^*T(\sigma)^{-1}T(z)x=0\) is solved at each step of the nonlinear iteration. When this flag is set the simpler form \(x^*T(z)x=0\) is used, which is supposed to be valid only for Hermitian problems.

Source code at slepc4py/SLEPc/NEP.pyx:2002

setRIIKSP(ksp)#

Set a linear solver object associated to the nonlinear eigensolver.

Collective.

Parameters:

ksp (KSP) – The linear solver object.

Return type:

None

Source code at slepc4py/SLEPc/NEP.pyx:2100

setRIILagPreconditioner(lag)#

Set when the preconditioner is rebuilt in the nonlinear solve.

Logically collective.

Parameters:

lag (int) – 0 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within the nonlinear iteration, 2 means every second time the Jacobian is built, etc.

Return type:

None

Source code at slepc4py/SLEPc/NEP.pyx:1878

setRIIMaximumIterations(its)#

Set the max. number of inner iterations to be used in the RII solver.

Logically collective.

These are the Newton iterations related to the computation of the nonlinear Rayleigh functional.

Parameters:

its (int) – Maximum inner iterations.

Return type:

None

Source code at slepc4py/SLEPc/NEP.pyx:1962

setRefine(ref, npart=None, tol=None, its=None, scheme=None)#

Set the refinement strategy used by the NEP object.

Logically collective.

Set the refinement strategy used by the NEP 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/NEP.pyx:694

setSLPDeflationThreshold(deftol)#

Set the threshold used to switch between deflated and non-deflated.

Logically collective.

Parameters:

deftol (float) – The threshold value.

Return type:

None

Notes

Normally, the solver iterates on the extended problem in order to deflate previously converged eigenpairs. If this threshold is set to a nonzero value, then once the residual error is below this threshold the solver will continue the iteration without deflation. The intention is to be able to improve the current eigenpair further, despite having previous eigenpairs with somewhat bad precision.

Source code at slepc4py/SLEPc/NEP.pyx:2139

setSLPEPS(eps)#

Set a linear eigensolver object associated to the nonlinear eigensolver.

Collective.

Parameters:

eps (EPS) – The linear eigensolver.

Return type:

None

Source code at slepc4py/SLEPc/NEP.pyx:2186

setSLPEPSLeft(eps)#

Set a linear eigensolver object associated to the nonlinear eigensolver.

Collective.

Used to compute left eigenvectors in the two-sided variant of SLP.

Parameters:

eps (EPS) – The linear eigensolver.

Return type:

None

Source code at slepc4py/SLEPc/NEP.pyx:2223

setSLPKSP(ksp)#

Set a linear solver object associated to the nonlinear eigensolver.

Collective.

Parameters:

ksp (KSP) – The linear solver object.

Return type:

None

Source code at slepc4py/SLEPc/NEP.pyx:2262

setSplitOperator(A, f, structure=None)#

Set the operator of the nonlinear eigenvalue problem in split form.

Collective.

Parameters:
Return type:

None

Notes

The nonlinear operator is written as \(T(\lambda) = \sum_i A_i f_i(\lambda)\), for \(i=1,\dots,n\). The derivative \(T'(\lambda)\) can be obtained using the derivatives of \(f_i\).

The structure flag provides information about \(A_i\)’s nonzero pattern.

This function must be called before setUp(). If it is called again after setUp() then the NEP object is reset.

Source code at slepc4py/SLEPc/NEP.pyx:1632

setSplitPreconditioner(P, structure=None)#

Set the operator in split form.

Collective.

Set the operator in split form from which to build the preconditioner to be used when solving the nonlinear eigenvalue problem in split form.

Parameters:
Return type:

None

Source code at slepc4py/SLEPc/NEP.pyx:1719

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/NEP.pyx:1029

Parameters:
Return type:

None

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/NEP.pyx:544

setTolerances(tol=None, maxit=None)#

Set the tolerance and max. iteration count used in convergence tests.

Logically collective.

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

  • maxit (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/NEP.pyx:596

setTrackAll(trackall)#

Set if the solver must compute the residual of all approximate eigenpairs.

Logically collective.

Parameters:

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

Return type:

None

Source code at slepc4py/SLEPc/NEP.pyx:777

setTwoSided(twosided)#

Set the solver to use a two-sided variant.

Logically collective.

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

Return type:

None

Notes

If the user sets twosided to True then the solver uses a variant of the algorithm that computes both right and left eigenvectors. This is usually much more costly. This option is not available in all solvers.

When using two-sided solvers, the problem matrices must have both the Mat.mult and Mat.multTranspose operations defined.

Source code at slepc4py/SLEPc/NEP.pyx:1803

setType(nep_type)#

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

Logically collective.

Parameters:

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

Return type:

None

Notes

The default is RII. Normally, it is best to use setFromOptions() and then set the NEP 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, NEPSetType

Source code at slepc4py/SLEPc/NEP.pyx:289

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 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, NEPSetUp

Source code at slepc4py/SLEPc/NEP.pyx:1174

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 NEP account for all the possible values. Also, some values make sense only for certain types of problems. If SLEPc is compiled for real numbers NEP.Which.LARGEST_IMAGINARY and NEP.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 NEP.Which.TARGET_IMAGINARY is available only in case PETSc and SLEPc have been built with complex scalars.

NEP.Which.ALL is intended for use in the context of the PEP.Type.CISS solver for computing all eigenvalues in a region.

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

solve()#

Solve the nonlinear eigenproblem.

Collective.

Notes

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/NEP.pyx:1195

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/NEP.pyx:1460

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/NEP.pyx:1479

view(viewer=None)#

Print the NEP data structure.

Collective.

Parameters:

viewer (Viewer | None) – Visualization context; if not provided, the standard output is used.

Return type:

None

See also

NEPView

Source code at slepc4py/SLEPc/NEP.pyx:223

Attributes Documentation

bv#

The basis vectors (BV) object associated.

Source code at slepc4py/SLEPc/NEP.pyx:3058

ds#

The direct solver (DS) object associated.

Source code at slepc4py/SLEPc/NEP.pyx:3072

max_it#

The maximum iteration count used by the NEP convergence tests.

Source code at slepc4py/SLEPc/NEP.pyx:3044

problem_type#

The problem type from the NEP object.

Source code at slepc4py/SLEPc/NEP.pyx:3016

rg#

The region (RG) object associated.

Source code at slepc4py/SLEPc/NEP.pyx:3065

target#

The value of the target.

Source code at slepc4py/SLEPc/NEP.pyx:3030

tol#

The tolerance used by the NEP convergence tests.

Source code at slepc4py/SLEPc/NEP.pyx:3037

track_all#

Compute the residual of all approximate eigenpairs.

Source code at slepc4py/SLEPc/NEP.pyx:3051

which#

The portion of the spectrum to be sought.

Source code at slepc4py/SLEPc/NEP.pyx:3023