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

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

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 necessary internal data structures.

setWhichEigenpairs(which)

Set which portion of the spectrum is to be sought.

solve()

Solve the eigensystem.

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:299

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

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

cancelMonitor()#

Clear all monitors for a NEP object.

Logically collective.

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

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

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

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

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

destroy()#

Destroy the NEP object.

Collective.

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

Return type:

Self

errorView(etype=None, viewer=None)#

Display the errors associated with the computed solution.

Collective.

Display the eigenvalues and the errors associated with the computed solution

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:1017

getBV()#

Get the basis vectors object associated to the eigensolver.

Not collective.

Returns:

The basis vectors context.

Return type:

BV

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

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:2078

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:2247

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:2229

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:2143

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:2192

getConverged()#

Get the number of converged eigenpairs.

Not collective.

Returns:

Number of converged eigenpairs.

Return type:

int

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

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:879

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:462

getDS()#

Get the direct solver associated to the eigensolver.

Not collective.

Returns:

The direct solver context.

Return type:

DS

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

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:608

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

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

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

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

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.

Parameters:
  • F – Function matrix

  • P – preconditioner matrix (usually the same as the F)

  • function – Function evaluation routine

Return type:

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

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

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:1809

getInterpolPEP()#

Get the associated polynomial eigensolver object.

Collective.

Get the polynomial eigensolver object associated with the nonlinear eigensolver.

Returns:

The polynomial eigensolver.

Return type:

PEP

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

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:861

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.

Parameters:
  • J – Jacobian matrix

  • jacobian – Jacobian evaluation routine

Return type:

tuple[petsc4py.PETSc.Mat, NEPJacobian]

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

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:939

getMonitor()#

Get the list of monitor functions.

Not collective.

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

Return type:

NEPMonitorFunction

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:1701

getNArnoldiLagPreconditioner()#

Get how often the preconditioner is rebuilt.

Not collective.

Returns:

The lag parameter.

Return type:

int

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

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:1988

getNLEIGSFullBasis()#

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

Not collective.

Returns:

True if the full-basis variant must be selected.

Return type:

bool

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

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:1922

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:2041

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:1885

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:2020

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:1849

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:269

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:325

getRG()#

Get the region object associated to the eigensolver.

Not collective.

Returns:

The region context.

Return type:

RG

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

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:1423

getRIIDeflationThreshold()#

Get the threshold value that controls deflation.

Not collective.

Returns:

The threshold value.

Return type:

float

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

getRIIHermitian()#

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

Not collective.

Get the flag about using the Hermitian version of the scalar nonlinear equation.

Returns:

If True, the Hermitian version is used.

Return type:

bool

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

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:1550

getRIILagPreconditioner()#

Get how often the preconditioner is rebuilt.

Not collective.

Returns:

The lag parameter.

Return type:

int

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

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:1455

getRefine()#

Get the refinement strategy used by the NEP object.

Not collective.

Get the refinement strategy used by the NEP object and the associated parameters.

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:493

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:563

getSLPDeflationThreshold()#

Get the threshold value that controls deflation.

Not collective.

Returns:

The threshold value.

Return type:

float

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

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:1610

getSLPEPSLeft()#

Get the left eigensolver.

Collective.

Returns:

The linear eigensolver.

Return type:

EPS

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

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:1670

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:1228

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:1289

getStoppingTest()#

Get the stopping function.

Not collective.

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

Return type:

NEPStoppingFunction

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:383

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:422

getTrackAll()#

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

Not collective.

Returns:

Whether the solver compute all residuals or not.

Return type:

bool

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

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:1314

getType()#

Get the NEP type of this object.

Not collective.

Returns:

The solver currently being used.

Return type:

str

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

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:354

reset()#

Reset the NEP object.

Collective.

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

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

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

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:2064

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:2210

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:2093

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:2173

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:478

setDS(ds)#

Set a direct solver object associated to the eigensolver.

Collective.

Parameters:

ds (DS) – The direct solver context.

Return type:

None

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

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

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

setFromOptions()#

Set NEP options from the options database.

Collective.

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:314

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:1078

setInitialSpace(space)#

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

Collective.

Parameters:

space (Vec) – The initial space

Return type:

None

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

setInterpolInterpolation(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 polynomial.

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:1787

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:1755

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:1140

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

Append a monitor function to the list of monitors.

Logically collective.

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

Parameters:
Return type:

None

setNArnoldiKSP(ksp)#

Set a linear solver object associated to the nonlinear eigensolver.

Collective.

Parameters:

ksp (petsc4py.PETSc.KSP) – The linear solver object.

Return type:

None

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

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:1717

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:1975

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

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

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:1900

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:1864

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

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

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:1829

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

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

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

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

setRG(rg)#

Set a region object associated to the eigensolver.

Collective.

Parameters:

rg (RG) – The region context.

Return type:

None

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

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

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

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

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

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

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

setRIIKSP(ksp)#

Set a linear solver object associated to the nonlinear eigensolver.

Collective.

Parameters:

ksp (petsc4py.PETSc.KSP) – The linear solver object.

Return type:

None

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

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:1378

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:1438

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 linear system solves

Return type:

None

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

setSLPDeflationThreshold(deftol)#

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

Logically collective.

Parameters:

deftol (float) – The threshold value.

Return type:

None

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

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:1597

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:1626

setSLPKSP(ksp)#

Set a linear solver object associated to the nonlinear eigensolver.

Collective.

Parameters:

ksp (petsc4py.PETSc.KSP) – The linear solver object.

Return type:

None

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

setSplitOperator(A, f, structure=None)#

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

Collective.

Parameters:
Return type:

None

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

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:1260

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:772

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

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

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

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

setTrackAll(trackall)#

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

Logically collective.

Parameters:

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

Return type:

None

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

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

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

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

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

setUp()#

Set up all the necessary internal data structures.

Collective.

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

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

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

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

solve()#

Solve the eigensystem.

Collective.

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

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:1046

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:1061

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

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

Attributes Documentation

bv#

The basis vectors (BV) object associated.

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

ds#

The direct solver (DS) object associated.

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

max_it#

The maximum iteration count used by the NEP convergence tests.

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

problem_type#

The problem type from the NEP object.

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

rg#

The region (RG) object associated.

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

target#

The value of the target.

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

tol#

The tolerance used by the NEP convergence tests.

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

track_all#

Compute the residual of all approximate eigenpairs.

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

which#

The portion of the spectrum to be sought.

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