slepc4py.SLEPc.NEP#
- class slepc4py.SLEPc.NEP#
Bases:
Object
NEP.
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
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.
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 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
NEP
object.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
KSP
object 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 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.
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 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
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:
- cancelMonitor()#
Clear all monitors for a
NEP
object.Logically collective.
Source code at slepc4py/SLEPc/NEP.pyx:821
- 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:
- create(comm=None)#
Create the NEP object.
Collective.
- destroy()#
Destroy the NEP object.
Collective.
Source code at slepc4py/SLEPc/NEP.pyx:194
- Return type:
- 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:
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_DETAIL
then a table with eigenvalues and corresponding errors is printed.
- 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:
- 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.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.
- getCISSRefinement()#
Get the values of various refinement parameters in the CISS solver.
Not collective.
- Returns:
- Return type:
- getCISSSizes()#
Get the values of various size parameters in the CISS solver.
Not collective.
- Returns:
- Return type:
- getCISSThreshold()#
Get the values of various threshold parameters in the CISS solver.
Not collective.
- Returns:
- Return type:
- getConverged()#
Get the number of converged eigenpairs.
Not collective.
- Returns:
Number of converged eigenpairs.
- Return type:
- getConvergedReason()#
Get the reason why the
solve()
iteration was stopped.Not collective.
- Returns:
Negative value indicates diverged, positive value converged.
- Return type:
- 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:
- 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:
- 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:
- 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:
- 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:
- getInterpolInterpolation()#
Get the tolerance and maximum degree for the interpolation polynomial.
Not collective.
- Returns:
- Return type:
- getInterpolPEP()#
Get the associated polynomial eigensolver object.
Collective.
Get the polynomial eigensolver object associated with the nonlinear eigensolver.
- Returns:
The polynomial eigensolver.
- Return type:
- 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:
- 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:
- getLeftEigenvector(i, Wr, Wi=None)#
Get the i-th left eigenvector as computed by
solve()
.Collective.
- Parameters:
- Return type:
Notes
The index
i
should be a value between0
andnconv-1
(seegetConverged()
). Eigensolutions are indexed according to the ordering criterion established withsetWhichEigenpairs()
.Left eigenvectors are available only if the twosided flag was set with
setTwoSided()
.
- getMonitor()#
Get the list of monitor functions.
Not collective.
Source code at slepc4py/SLEPc/NEP.pyx:813
- Return type:
- getNArnoldiKSP()#
Get the linear solver object associated with the nonlinear eigensolver.
Collective.
- Returns:
The linear solver object.
- Return type:
- 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:
- 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:
- 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.KSP
solvers is equal to the number of shifts provided by the user, or 1 if the user did not provide shifts.
- getNLEIGSLocking()#
Get the locking flag used in the NLEIGS method.
Not collective.
- Returns:
The locking flag.
- Return type:
- getNLEIGSRKShifts()#
Get the list of shifts used in the Rational Krylov method.
Not collective.
- Returns:
The shift values.
- Return type:
- getNLEIGSRestart()#
Get the restart parameter used in the NLEIGS method.
Not collective.
- Returns:
The number of vectors to be kept at restart.
- Return type:
- 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
NEP
object.Not collective.
- Returns:
The problem type that was previously set.
- Return type:
- 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, the
petsc4py.PETSc.KSP
relative 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.
Get the flag about using the Hermitian version of the scalar nonlinear equation.
- Returns:
If True, the Hermitian version is used.
- Return type:
- getRIIKSP()#
Get the linear solver object associated with the nonlinear eigensolver.
Collective.
- Returns:
The linear solver object.
- Return type:
- 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.
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:
- getRefineKSP()#
Get the
KSP
object used by the eigensolver in the refinement phase.Collective.
- Returns:
The linear solver object.
- Return type:
- 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:
- getSLPEPSLeft()#
Get the left eigensolver.
Collective.
- Returns:
The linear eigensolver.
- Return type:
- getSLPKSP()#
Get the linear solver object associated with the nonlinear eigensolver.
Collective.
- Returns:
The linear solver object.
- Return type:
- getSplitOperator()#
Get the operator of the nonlinear eigenvalue problem in split form.
Collective.
- Returns:
A (
list
ofpetsc4py.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 (
list
ofpetsc4py.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 function.
Not collective.
Source code at slepc4py/SLEPc/NEP.pyx:782
- Return type:
- 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.
- 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:
- 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:
- 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:
- getType()#
Get the NEP type of this object.
Not collective.
- Returns:
The solver currently being used.
- Return type:
- 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:
- reset()#
Reset the NEP object.
Collective.
Source code at slepc4py/SLEPc/NEP.pyx:204
- 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:
- setCISSRefinement(inner=None, blsize=None)#
Set the values of various refinement parameters in the CISS solver.
Logically collective.
- Parameters:
- Return type:
- 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.KSP
object is shared among all processes of theNEP
communicator. Otherwise, the communicator is split into npart communicators, so thatnpart
petsc4py.PETSc.KSP
solves proceed simultaneously.
- setCISSThreshold(delta=None, spur=None)#
Set the values of various threshold parameters in the CISS solver.
Logically collective.
- Parameters:
- Return type:
- 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:
- 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:
- 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:304
- 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:
- setInitialSpace(space)#
Set the initial space from which the eigensolver starts to iterate.
Collective.
- 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:
- Return type:
- setInterpolPEP(pep)#
Set a polynomial eigensolver object associated to the nonlinear eigensolver.
Collective.
- 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:
- setMonitor(monitor, args=None, kargs=None)#
Append a monitor function to the list of monitors.
Logically collective.
- setNArnoldiKSP(ksp)#
Set a linear solver object associated to the nonlinear eigensolver.
Collective.
- Parameters:
ksp (petsc4py.PETSc.KSP) – The linear solver object.
- Return type:
- 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.
- 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.
- 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 behaviour can be changed so that all directions are kept in the working subspace even if already converged to working accuracy (the non-locking variant).
- setNLEIGSRKShifts(shifts)#
Set a list of shifts to be used in the Rational Krylov method.
Collective.
- 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.
- setOptionsPrefix(prefix=None)#
Set the prefix used for searching for all NEP options in the database.
Logically collective.
- setProblemType(problem_type)#
Set the type of the eigenvalue problem.
Logically collective.
- Parameters:
problem_type (ProblemType) – The problem type to be set.
- Return type:
- 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, the
petsc4py.PETSc.KSP
relative tolerance is constant.- Return type:
- 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.
- 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.
- setRIIKSP(ksp)#
Set a linear solver object associated to the nonlinear eigensolver.
Collective.
- Parameters:
ksp (petsc4py.PETSc.KSP) – The linear solver object.
- Return type:
- 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:
- setSLPDeflationThreshold(deftol)#
Set the threshold used to switch between deflated and non-deflated.
Logically collective.
- setSLPEPS(eps)#
Set a linear eigensolver object associated to the nonlinear eigensolver.
Collective.
- 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.
- setSLPKSP(ksp)#
Set a linear solver object associated to the nonlinear eigensolver.
Collective.
- Parameters:
ksp (petsc4py.PETSc.KSP) – The linear solver object.
- Return type:
- 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:
- 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.
- 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()
.
- setTolerances(tol=None, maxit=None)#
Set the tolerance and max. iteration count used in convergence tests.
Logically collective.
- Parameters:
- Return type:
- setTrackAll(trackall)#
Set if the solver must compute the residual of all approximate eigenpairs.
Logically collective.
- 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.
- setType(nep_type)#
Set the particular solver to be used in the NEP object.
Logically collective.
- 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:832
- 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:
- solve()#
Solve the eigensystem.
Collective.
Source code at slepc4py/SLEPc/NEP.pyx:843
- Return type:
- valuesView(viewer=None)#
Display the computed eigenvalues in a viewer.
Collective.
- vectorsView(viewer=None)#
Output computed eigenvectors to a viewer.
Collective.
- view(viewer=None)#
Print the NEP data structure.
Collective.
Attributes Documentation
- 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.