slepc4py.SLEPc.NEP#
- class slepc4py.SLEPc.NEP#
Bases:
ObjectNonlinear 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
NEP CISS extraction technique.
NEP convergence test.
NEP convergence reasons.
NEP error type to assess accuracy of computed solutions.
NEP problem type.
NEP refinement strategy.
NEP scheme for solving linear systems during iterative refinement.
NEP stopping test.
NEP type.
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.
Clear all monitors for a
NEPobject.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.
Get the extraction technique used in the CISS solver.
Get the list of linear solver objects associated with the CISS solver.
Get the values of various refinement parameters in the CISS solver.
Get the values of various size parameters in the CISS solver.
Get the values of various threshold parameters in the CISS solver.
Get the number of converged eigenpairs.
Get the reason why the
solve()iteration was stopped.Get the method used to compute the error estimate used in the convergence test.
getDS()Get the direct solver associated to the eigensolver.
Get the number of eigenvalues to compute.
getEigenpair(i[, Vr, Vi])Get the i-th solution of the eigenproblem as computed by
solve().Get the eigenvalue comparison function.
Get the error estimate associated to the i-th computed eigenpair.
Get the function to compute the nonlinear Function \(T(\lambda)\).
Get the tolerance and maximum degree for the interpolation polynomial.
Get the associated polynomial eigensolver object.
Get the current iteration number.
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().Get the list of monitor functions.
Get the linear solver object associated with the nonlinear eigensolver.
Get how often the preconditioner is rebuilt.
Get the linear eigensolver object associated with the nonlinear eigensolver.
Get the flag that indicates if NLEIGS is using the full-basis variant.
Get the tolerance and maximum degree for the interpolation polynomial.
Get the list of linear solver objects associated with the NLEIGS solver.
Get the locking flag used in the NLEIGS method.
Get the list of shifts used in the Rational Krylov method.
Get the restart parameter used in the NLEIGS method.
Get the prefix used for searching for all NEP options in the database.
Get the problem type from the
NEPobject.getRG()Get the region object associated to the eigensolver.
Get the constant tolerance flag.
Get the threshold value that controls deflation.
Get if the Hermitian version must be used by the solver.
Get the linear solver object associated with the nonlinear eigensolver.
Get how often the preconditioner is rebuilt.
Get the maximum number of inner iterations of RII.
Get the refinement strategy used by the NEP object.
Get the
KSPobject used by the eigensolver in the refinement phase.Get the threshold value that controls deflation.
Get the linear eigensolver object associated with the nonlinear eigensolver.
Get the left eigensolver.
Get the linear solver object associated with the nonlinear eigensolver.
Get the operator of the nonlinear eigenvalue problem in split form.
Get the operator of the split preconditioner.
Get the stopping test function.
Get the value of the target.
Get the tolerance and maximum iteration count.
Get the flag indicating whether all residual norms must be computed.
Get the flag indicating if a two-sided variant is being used.
getType()Get the NEP type of this object.
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.
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.
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.
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.
Set when the preconditioner is rebuilt in the nonlinear solve.
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
The basis vectors (
BV) object associated.The direct solver (
DS) object associated.The maximum iteration count used by the NEP convergence tests.
The problem type from the NEP object.
The region (
RG) object associated.The value of the target.
The tolerance used by the NEP convergence tests.
Compute the residual of all approximate eigenpairs.
The portion of the spectrum to be sought.
Methods Documentation
- appendOptionsPrefix(prefix=None)#
Append to the prefix used for searching for all NEP options in the database.
Logically collective.
- applyResolvent(omega, v, r, rg=None)#
Apply the resolvent \(T^{-1}(z)\) to a given vector.
Collective.
- Parameters:
- Return type:
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 regionrgis given then only those corresponding to eigenvalues inside the region are considered.See also
- cancelMonitor()#
Clear all monitors for a
NEPobject.Logically collective.
See also
Source code at slepc4py/SLEPc/NEP.pyx:1159
- Return type:
- computeError(i, etype=None)#
Compute the error associated with the i-th computed eigenpair.
Collective.
Compute the error (based on the residual norm) associated with the i-th computed eigenpair.
- Parameters:
- Returns:
The error bound, computed in various ways from the residual norm \(\|T(\lambda)x\|_2\) where \(\lambda\) is the eigenvalue and \(x\) is the eigenvector.
- Return type:
Notes
The index
ishould be a value between0andnconv-1(seegetConverged()).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.See also
- create(comm=None)#
Create the NEP object.
Collective.
- Parameters:
comm (Comm | None) – MPI communicator. If not provided, it defaults to all processes.
- Return type:
See also
- destroy()#
Destroy the NEP object.
Collective.
See also
Source code at slepc4py/SLEPc/NEP.pyx:242
- Return type:
- errorView(etype=None, viewer=None)#
Display the errors associated with the computed solution.
Collective.
Display the errors and the eigenvalues.
- Parameters:
viewer (petsc4py.PETSc.Viewer | None) – Visualization context; if not provided, the standard output is used.
- Return type:
Notes
By default, this function checks the error of all eigenpairs and prints the eigenvalues if all of them are below the requested tolerance. If the viewer has format
ASCII_INFO_DETAILthen a table with eigenvalues and corresponding errors is printed.See also
- getBV()#
Get the basis vectors object associated to the eigensolver.
Not collective.
- Returns:
The basis vectors context.
- Return type:
- getCISSExtraction()#
Get the extraction technique used in the CISS solver.
Not collective.
- Returns:
The extraction technique.
- Return type:
See also
- getCISSKSPs()#
Get the list of linear solver objects associated with the CISS solver.
Collective.
- Returns:
The linear solver objects.
- Return type:
Notes
The number of
petsc4py.PETSc.KSPsolvers is equal to the number of integration points divided by the number of partitions. This value is halved in the case of real matrices with a region centered at the real axis.See also
- getCISSRefinement()#
Get the values of various refinement parameters in the CISS solver.
Not collective.
- Returns:
- Return type:
See also
- getCISSSizes()#
Get the values of various size parameters in the CISS solver.
Not collective.
- Returns:
- Return type:
See also
- getCISSThreshold()#
Get the values of various threshold parameters in the CISS solver.
Not collective.
- Returns:
- Return type:
See also
- getConverged()#
Get the number of converged eigenpairs.
Not collective.
- Returns:
nconv – Number of converged eigenpairs.
- Return type:
Notes
This function should be called after
solve()has finished.The value
nconvmay be different from the number of requested solutionsnev, but not larger thanncv, seesetDimensions().See also
- getConvergedReason()#
Get the reason why the
solve()iteration was stopped.Not collective.
- Returns:
Negative value indicates diverged, positive value converged.
- Return type:
See also
- getConvergenceTest()#
Get the method used to compute the error estimate used in the convergence test.
Not collective.
- Returns:
The method used to compute the error estimate used in the convergence test.
- Return type:
See also
- getDS()#
Get the direct solver associated to the eigensolver.
Not collective.
- Returns:
The direct solver context.
- Return type:
- getDimensions()#
Get the number of eigenvalues to compute.
Not collective.
Get the number of eigenvalues to compute, and the dimension of the subspace.
- Returns:
- Return type:
See also
- getEigenpair(i, Vr=None, Vi=None)#
Get the i-th solution of the eigenproblem as computed by
solve().Collective.
The solution consists of both the eigenvalue and the eigenvector.
- Parameters:
- Returns:
The computed eigenvalue.
- Return type:
Notes
The index
ishould be a value between0andnconv-1(seegetConverged()). Eigenpairs are indexed according to the ordering criterion established withsetWhichEigenpairs().The eigenvector is normalized to have unit norm.
See also
- getEigenvalueComparison()#
Get the eigenvalue comparison function.
Not collective.
- Returns:
The eigenvalue comparison function.
- Return type:
See also
- getErrorEstimate(i)#
Get the error estimate associated to the i-th computed eigenpair.
Not collective.
- Parameters:
i (int) – Index of the solution to be considered.
- Returns:
Error estimate.
- Return type:
Notes
This is the error estimate used internally by the eigensolver. The actual error bound can be computed with
computeError().See also
- 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:
F (
petsc4py.PETSc.Mat) – Function matrix.P (
petsc4py.PETSc.Mat) – Preconditioner matrix (usually the same as the F).function (
NEPFunction) – Function evaluation routine.
- Return type:
See also
- getInterpolInterpolation()#
Get the tolerance and maximum degree for the interpolation polynomial.
Not collective.
- Returns:
- Return type:
- getInterpolPEP()#
Get the associated polynomial eigensolver object.
Collective.
- Returns:
The polynomial eigensolver.
- Return type:
See also
- getIterationNumber()#
Get the current iteration number.
Not collective.
If the call to
solve()is complete, then it returns the number of iterations carried out by the solution method.- Returns:
Iteration number.
- Return type:
See also
- 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:
J (
petsc4py.PETSc.Mat) – Jacobian matrix.jacobian (
NEPJacobian) – Jacobian evaluation routine.
- Return type:
See also
- getLeftEigenvector(i, Wr, Wi=None)#
Get the i-th left eigenvector as computed by
solve().Collective.
- Parameters:
- Return type:
Notes
The index
ishould be a value between0andnconv-1(seegetConverged()). Eigensolutions are indexed according to the ordering criterion established withsetWhichEigenpairs().Left eigenvectors are available only if the
twosidedflag was set withsetTwoSided().See also
getEigenpair,getConverged,setTwoSided,NEPGetLeftEigenvector
- getMonitor()#
Get the list of monitor functions.
Not collective.
- Returns:
The list of monitor functions.
- Return type:
See also
- getNArnoldiKSP()#
Get the linear solver object associated with the nonlinear eigensolver.
Collective.
- Returns:
The linear solver object.
- Return type:
See also
- getNArnoldiLagPreconditioner()#
Get how often the preconditioner is rebuilt.
Not collective.
- Returns:
The lag parameter.
- Return type:
- getNLEIGSEPS()#
Get the linear eigensolver object associated with the nonlinear eigensolver.
Collective.
- Returns:
The linear eigensolver.
- Return type:
See also
- getNLEIGSFullBasis()#
Get the flag that indicates if NLEIGS is using the full-basis variant.
Not collective.
- Returns:
Trueif the full-basis variant is selected.- Return type:
See also
- 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:
- Return type:
- getNLEIGSKSPs()#
Get the list of linear solver objects associated with the NLEIGS solver.
Collective.
- Returns:
The linear solver objects.
- Return type:
Notes
The number of
petsc4py.PETSc.KSPsolvers is equal to the number of shifts provided by the user, or 1 if the user did not provide shifts.See also
- getNLEIGSLocking()#
Get the locking flag used in the NLEIGS method.
Not collective.
- Returns:
The locking flag.
- Return type:
See also
- getNLEIGSRKShifts()#
Get the list of shifts used in the Rational Krylov method.
Not collective.
- Returns:
The shift values.
- Return type:
See also
- getNLEIGSRestart()#
Get the restart parameter used in the NLEIGS method.
Not collective.
- Returns:
The number of vectors to be kept at restart.
- Return type:
See also
- 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:
- getProblemType()#
Get the problem type from the
NEPobject.Not collective.
- Returns:
The problem type that was previously set.
- Return type:
See also
- getRG()#
Get the region object associated to the eigensolver.
Not collective.
- Returns:
The region context.
- Return type:
- getRIIConstCorrectionTol()#
Get the constant tolerance flag.
Not collective.
- Returns:
If
True, thepetsc4py.PETSc.KSPrelative tolerance is constant.- Return type:
- getRIIDeflationThreshold()#
Get the threshold value that controls deflation.
Not collective.
- Returns:
The threshold value.
- Return type:
- getRIIHermitian()#
Get if the Hermitian version must be used by the solver.
Not collective.
- Returns:
If
True, the Hermitian version is used.- Return type:
See also
- getRIIKSP()#
Get the linear solver object associated with the nonlinear eigensolver.
Collective.
- Returns:
The linear solver object.
- Return type:
See also
- getRIILagPreconditioner()#
Get how often the preconditioner is rebuilt.
Not collective.
- Returns:
The lag parameter.
- Return type:
- getRIIMaximumIterations()#
Get the maximum number of inner iterations of RII.
Not collective.
- Returns:
Maximum inner iterations.
- Return type:
- 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:
See also
- getRefineKSP()#
Get the
KSPobject used by the eigensolver in the refinement phase.Collective.
- Returns:
The linear solver object.
- Return type:
See also
- getSLPDeflationThreshold()#
Get the threshold value that controls deflation.
Not collective.
- Returns:
The threshold value.
- Return type:
- getSLPEPS()#
Get the linear eigensolver object associated with the nonlinear eigensolver.
Collective.
- Returns:
The linear eigensolver.
- Return type:
See also
- getSLPEPSLeft()#
Get the left eigensolver.
Collective.
- Returns:
The linear eigensolver.
- Return type:
See also
- getSLPKSP()#
Get the linear solver object associated with the nonlinear eigensolver.
Collective.
- Returns:
The linear solver object.
- Return type:
See also
- getSplitOperator()#
Get the operator of the nonlinear eigenvalue problem in split form.
Collective.
- Returns:
A (
listofpetsc4py.PETSc.Mat) – Coefficient matrices of the split form.structure (
petsc4py.PETSc.Mat.Structure) – Structure flag for matrices.
- Return type:
tuple[list[petsc4py.PETSc.Mat], list[FN], petsc4py.PETSc.Mat.Structure]
- getSplitPreconditioner()#
Get the operator of the split preconditioner.
Not collective.
- Returns:
P (
listofpetsc4py.PETSc.Mat) – Coefficient matrices of the split preconditioner.structure (
petsc4py.PETSc.Mat.Structure) – Structure flag for matrices.
- Return type:
tuple[list[petsc4py.PETSc.Mat], petsc4py.PETSc.Mat.Structure]
- getStoppingTest()#
Get the stopping test function.
Not collective.
- Returns:
The stopping test function.
- Return type:
See also
- getTarget()#
Get the value of the target.
Not collective.
- Returns:
The value of the target.
- Return type:
Notes
If the target was not set by the user, then zero is returned.
See also
- getTolerances()#
Get the tolerance and maximum iteration count.
Not collective.
Get the tolerance and maximum iteration count used by the default NEP convergence tests.
- Returns:
- Return type:
See also
- getTrackAll()#
Get the flag indicating whether all residual norms must be computed.
Not collective.
- Returns:
Whether the solver computes all residuals or not.
- Return type:
See also
- 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:
See also
- getType()#
Get the NEP type of this object.
Not collective.
- Returns:
The solver currently being used.
- Return type:
See also
- getWhichEigenpairs()#
Get which portion of the spectrum is to be sought.
Not collective.
- Returns:
The portion of the spectrum to be sought by the solver.
- Return type:
See also
- reset()#
Reset the NEP object.
Collective.
See also
Source code at slepc4py/SLEPc/NEP.pyx:256
- Return type:
- setBV(bv)#
Set the basis vectors object associated to the eigensolver.
Collective.
- setCISSExtraction(extraction)#
Set the extraction technique used in the CISS solver.
Logically collective.
- Parameters:
extraction (CISSExtraction) – The extraction technique.
- Return type:
See also
- setCISSRefinement(inner=None, blsize=None)#
Set the values of various refinement parameters in the CISS solver.
Logically collective.
- Parameters:
- Return type:
See also
- setCISSSizes(ip=None, bs=None, ms=None, npart=None, bsmax=None, realmats=False)#
Set the values of various size parameters in the CISS solver.
Logically collective.
- Parameters:
- Return type:
Notes
The default number of partitions is 1. This means the internal
petsc4py.PETSc.KSPobject is shared among all processes of theNEPcommunicator. Otherwise, the communicator is split intonpartcommunicators, so thatnpartpetsc4py.PETSc.KSPsolves proceed simultaneously.See also
getCISSSizes,setCISSThreshold,setCISSRefinement,NEPCISSSetSizes
- setCISSThreshold(delta=None, spur=None)#
Set the values of various threshold parameters in the CISS solver.
Logically collective.
- Parameters:
- Return type:
See also
- setConvergenceTest(conv)#
Set how to compute the error estimate used in the convergence test.
Logically collective.
- Parameters:
conv (Conv) – The method used to compute the error estimate used in the convergence test.
- Return type:
See also
- setDS(ds)#
Set a direct solver object associated to the eigensolver.
Collective.
- setDimensions(nev=None, ncv=None, mpd=None)#
Set the number of eigenvalues to compute.
Logically collective.
Set the number of eigenvalues to compute and the dimension of the subspace.
- Parameters:
- Return type:
Notes
Use
DETERMINEforncvandmpdto assign a reasonably good value, which is dependent on the solution method.The parameters
ncvandmpdare intimately related, so that the user is advised to set one of them at most. Normal usage is the following:In cases where
nevis small, the user setsncv(a reasonable default is 2 *nev).In cases where
nevis large, the user setsmpd.
The value of
ncvshould always be betweennevand (nev+mpd), typicallyncv=nev+mpd. Ifnevis not too large,mpd=nevis a reasonable choice, otherwise a smaller value should be used.See also
- setEigenvalueComparison(comparison, args=None, kargs=None)#
Set an eigenvalue comparison function.
Logically collective.
Notes
This eigenvalue comparison function is used when
setWhichEigenpairs()is set toNEP.Which.USER.
- setFromOptions()#
Set NEP options from the options database.
Collective.
Notes
To see all options, run your program with the
-helpoption.This routine must be called before
setUp()if the user is to be allowed to set the solver type.See also
Source code at slepc4py/SLEPc/NEP.pyx:404
- Return type:
- 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:
function (NEPFunction) – Function evaluation routine.
F (petsc4py.PETSc.Mat | None) – Function matrix.
P (petsc4py.PETSc.Mat | None) – Preconditioner matrix (usually the same as
F).
- Return type:
See also
- setInitialSpace(space)#
Set the initial space from which the eigensolver starts to iterate.
Collective.
Notes
Some solvers start to iterate on a single vector (initial vector). In that case, only the first vector is taken into account and the other vectors are ignored.
These vectors do not persist from one
solve()call to the other, so the initial space should be set every time.The vectors do not need to be mutually orthonormal, since they are explicitly orthonormalized internally.
Common usage of this function is when the user can provide a rough approximation of the wanted eigenspace. Then, convergence may be faster.
See also
- setInterpolInterpolation(tol=None, deg=None)#
Set the tolerance and maximum degree for the interpolation polynomial.
Collective.
- Parameters:
- Return type:
- setInterpolPEP(pep)#
Set a polynomial eigensolver object associated to the nonlinear eigensolver.
Collective.
See also
- 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:
jacobian (NEPJacobian) – Jacobian evaluation routine.
J (petsc4py.PETSc.Mat | None) – Jacobian matrix.
- Return type:
See also
- setMonitor(monitor, args=None, kargs=None)#
Append a monitor function to the list of monitors.
Logically collective.
See also
- setNArnoldiKSP(ksp)#
Set a linear solver object associated to the nonlinear eigensolver.
Collective.
See also
- 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:
Notes
The default is 1. The preconditioner is ALWAYS built in the first iteration of a nonlinear solve.
- setNLEIGSEPS(eps)#
Set a linear eigensolver object associated to the nonlinear eigensolver.
Collective.
See also
- 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.
Notes
The default is to use a compact representation of the Krylov basis, that is, \(V = (I \otimes U) S\), with a
BVof typeTENSOR. 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, seesetTwoSided().In the full-basis variant, the NLEIGS solver uses an
EPSobject to explicitly solve the linearized eigenproblem, seegetNLEIGSEPS().
- 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:
- Return type:
- setNLEIGSLocking(lock)#
Toggle between locking and non-locking variants of the NLEIGS method.
Logically collective.
Notes
The default is to lock converged eigenpairs when the method restarts. This behavior can be changed so that all directions are kept in the working subspace even if already converged to working accuracy (the non-locking variant).
See also
- setNLEIGSRKShifts(shifts)#
Set a list of shifts to be used in the Rational Krylov method.
Collective.
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.
See also
- setNLEIGSRestart(keep)#
Set the restart parameter for the NLEIGS method.
Logically collective.
The proportion of basis vectors that must be kept after restart.
Notes
Allowed values are in the range [0.1,0.9]. The default is 0.5.
See also
- 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:
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_")
- setProblemType(problem_type)#
Set the type of the eigenvalue problem.
Logically collective.
- Parameters:
problem_type (ProblemType) – The problem type to be set.
- Return type:
Notes
This function is used to provide a hint to the
NEPsolver 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.See also
- setRG(rg)#
Set a region object associated to the eigensolver.
Collective.
- setRIIConstCorrectionTol(cct)#
Set a flag to keep the tolerance used in the linear solver constant.
Logically collective.
- Parameters:
cct (bool) – If
True, thepetsc4py.PETSc.KSPrelative tolerance is constant.- Return type:
Notes
By default, an exponentially decreasing tolerance is set in the
KSPused within the nonlinear iteration, so that each Newton iteration requests better accuracy than the previous one. The constant correction tolerance flag stops this behavior.
- 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.
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.
- 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.
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.
See also
- setRIIKSP(ksp)#
Set a linear solver object associated to the nonlinear eigensolver.
Collective.
See also
- 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:
- 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.
- 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:
- Return type:
See also
- setSLPDeflationThreshold(deftol)#
Set the threshold used to switch between deflated and non-deflated.
Logically collective.
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.
- setSLPEPS(eps)#
Set a linear eigensolver object associated to the nonlinear eigensolver.
Collective.
See also
- 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.
See also
- setSLPKSP(ksp)#
Set a linear solver object associated to the nonlinear eigensolver.
Collective.
See also
- setSplitOperator(A, f, structure=None)#
Set the operator of the nonlinear eigenvalue problem in split form.
Collective.
- Parameters:
A (petsc4py.PETSc.Mat | list[petsc4py.PETSc.Mat]) – Coefficient matrices of the split form.
structure (petsc4py.PETSc.Mat.Structure | None) – Structure flag for matrices.
- Return type:
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
structureflag provides information about \(A_i\)’s nonzero pattern.This function must be called before
setUp(). If it is called again aftersetUp()then theNEPobject is reset.See also
- 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:
P (petsc4py.PETSc.Mat | list[petsc4py.PETSc.Mat]) – Coefficient matrices of the split preconditioner.
structure (petsc4py.PETSc.Mat.Structure | None) – Structure flag for matrices.
- Return type:
- setStoppingTest(stopping, args=None, kargs=None)#
Set a function to decide when to stop the outer iteration of the eigensolver.
Logically collective.
See also
- setTarget(target)#
Set the value of the target.
Logically collective.
Notes
The target is a scalar value used to determine the portion of the spectrum of interest. It is used in combination with
setWhichEigenpairs().When PETSc is built with real scalars, it is not possible to specify a complex target.
See also
- setTolerances(tol=None, maxit=None)#
Set the tolerance and max. iteration count used in convergence tests.
Logically collective.
- Parameters:
- Return type:
Notes
Use
DETERMINEformax_itto assign a reasonably good value, which is dependent on the solution method.See also
- setTrackAll(trackall)#
Set if the solver must compute the residual of all approximate eigenpairs.
Logically collective.
See also
- 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.
Notes
If the user sets
twosidedtoTruethen 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.multandMat.multTransposeoperations defined.See also
- setType(nep_type)#
Set the particular solver to be used in the NEP object.
Logically collective.
Notes
The default is
RII. Normally, it is best to usesetFromOptions()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
- 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.Source code at slepc4py/SLEPc/NEP.pyx:1174
- Return type:
- setWhichEigenpairs(which)#
Set which portion of the spectrum is to be sought.
Logically collective.
- Parameters:
which (Which) – The portion of the spectrum to be sought by the solver.
- Return type:
Notes
Not all eigensolvers implemented in 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_IMAGINARYandNEP.Which.SMALLEST_IMAGINARYuse the absolute value of the imaginary part for eigenvalue selection.The target is a scalar value provided with
setTarget().The criterion
NEP.Which.TARGET_IMAGINARYis available only in case PETSc and SLEPc have been built with complex scalars.NEP.Which.ALLis intended for use in the context of thePEP.Type.CISSsolver for computing all eigenvalues in a region.See also
- solve()#
Solve the nonlinear eigenproblem.
Collective.
Notes
solve()will return without generating an error regardless of whether all requested solutions were computed or not. CallgetConverged()to get the actual number of computed solutions, andgetConvergedReason()to determine if the solver converged or failed and why.See also
Source code at slepc4py/SLEPc/NEP.pyx:1195
- Return type:
- valuesView(viewer=None)#
Display the computed eigenvalues in a viewer.
Collective.
- Parameters:
viewer (Viewer | None) – Visualization context; if not provided, the standard output is used.
- Return type:
See also
- vectorsView(viewer=None)#
Output computed eigenvectors to a viewer.
Collective.
- Parameters:
viewer (Viewer | None) – Visualization context; if not provided, the standard output is used.
- Return type:
See also
- view(viewer=None)#
Print the NEP data structure.
Collective.
- Parameters:
viewer (Viewer | None) – Visualization context; if not provided, the standard output is used.
- Return type:
See also
Attributes Documentation
- max_it#
The maximum iteration count used by the NEP convergence tests.
- problem_type#
The problem type from the NEP object.
- 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.