slepc4py.SLEPc.SVD#
- class slepc4py.SLEPc.SVD#
Bases:
ObjectSingular Value Decomposition Solver.
The Singular Value Decomposition Solver (
SVD) is very similar to theEPSobject, but intended for the computation of the partial SVD of a rectangular matrix. With this type of object, the user can specify an SVD problem and solve it with any of the different solvers encapsulated by the package. Some of these solvers are actually implemented through calls toEPSeigensolvers.Enumerations
SVD convergence test.
SVD convergence reasons.
SVD error type to assess accuracy of computed solutions.
SVD problem type.
SVD stopping test.
SVD TRLanczos bidiagonalization choices for the GSVD case.
SVD type.
SVD desired part of spectrum.
Methods Summary
appendOptionsPrefix([prefix])Append to the prefix used for searching for all SVD options in the database.
Clear all monitors for an
SVDobject.computeError(i[, etype])Compute the error associated with the i-th singular triplet.
create([comm])Create the SVD object.
destroy()Destroy the SVD object.
errorView([etype, viewer])Display the errors associated with the computed solution.
getBV()Get the basis vectors objects associated to the SVD object.
Get the number of converged singular triplets.
Get the reason why the
solve()iteration was stopped.Get the method used to compute the error estimate used in the convergence test.
Get the eigensolver object associated to the singular value solver.
Get the flag indicating if \(A^*A\) is built explicitly.
Get the eigensolver object associated to the singular value solver.
Get the flag indicating if \(H(A)\) is built explicitly.
getDS()Get the direct solver associated to the singular value solver.
Get the number of singular values to compute and the dimension of the subspace.
Get the mode used to handle the transpose of the associated matrix.
Get the current iteration number.
Get if the variant of the Lanczos method to be used is one-sided or two-sided.
Get the list of monitor functions.
Get the matrices associated with the singular value problem.
Get the prefix used for searching for all SVD options in the database.
Get the problem type from the SVD object.
getSignature([omega])Get the signature matrix defining a hyperbolic singular value problem.
getSingularTriplet(i[, U, V])Get the i-th triplet of the singular value decomposition.
Get the stopping test function.
Get the flag indicating if \(Z=[A^*,B^*]^*\) is built explicitly.
Get bidiagonalization choice used in the GSVD TRLanczos solver.
Get the linear solver object associated with the SVD solver.
Get the locking flag used in the thick-restart Lanczos method.
Get if the variant of the method to be used is one-sided or two-sided.
Get the restart parameter used in the thick-restart Lanczos method.
Get the threshold used in the threshold stopping test.
Get the tolerance and maximum iteration count.
Get the flag indicating if all residual norms must be computed or not.
getType()Get the SVD type of this object.
getValue(i)Get the i-th singular value as computed by
solve().getVectors(i, U, V)Get the i-th left and right singular vectors as computed by
solve().Get which singular triplets are to be sought.
Tell if the SVD corresponds to a generalized singular value problem.
Tell whether the SVD object corresponds to a hyperbolic singular value problem.
reset()Reset the SVD object.
setBV(V[, U])Set basis vectors objects associated to the SVD solver.
setConvergenceTest(conv)Set how to compute the error estimate used in the convergence test.
setCrossEPS(eps)Set an eigensolver object associated to the singular value solver.
setCrossExplicitMatrix([flag])Set if the eigensolver operator \(A^*A\) must be computed.
setCyclicEPS(eps)Set an eigensolver object associated to the singular value solver.
setCyclicExplicitMatrix([flag])Set if the eigensolver operator \(H(A)\) must be computed explicitly.
setDS(ds)Set a direct solver object associated to the singular value solver.
setDimensions([nsv, ncv, mpd])Set the number of singular values to compute and the dimension of the subspace.
Set SVD options from the options database.
setImplicitTranspose(mode)Set how to handle the transpose of the associated matrix.
setInitialSpace([spaceright, spaceleft])Set the initial spaces from which the SVD solver starts to iterate.
setLanczosOneSide([flag])Set if the variant of the Lanczos method to be used is one-sided or two-sided.
setMonitor(monitor[, args, kargs])Append a monitor function to the list of monitors.
setOperators(A[, B])Set the matrices associated with the singular value problem.
setOptionsPrefix([prefix])Set the prefix used for searching for all SVD options in the database.
setProblemType(problem_type)Set the type of the singular value problem.
setSignature([omega])Set the signature matrix defining a hyperbolic singular value problem.
setStoppingTest(stopping[, args, kargs])Set a function to decide when to stop the outer iteration of the eigensolver.
setTRLanczosExplicitMatrix([flag])Set if the matrix \(Z=[A^*,B^*]^*\) must be built explicitly.
setTRLanczosGBidiag(bidiag)Set the bidiagonalization choice to use in the GSVD TRLanczos solver.
setTRLanczosKSP(ksp)Set a linear solver object associated to the SVD solver.
setTRLanczosLocking(lock)Toggle between locking and non-locking variants of TRLanczos.
setTRLanczosOneSide([flag])Set if the variant of the method to be used is one-sided or two-sided.
setTRLanczosRestart(keep)Set the restart parameter for the thick-restart Lanczos method.
setThreshold(thres[, rel])Set the threshold used in the threshold stopping test.
setTolerances([tol, max_it])Set the tolerance and maximum iteration count used.
setTrackAll(trackall)Set flag to compute the residual of all singular triplets.
setType(svd_type)Set the particular solver to be used in the SVD object.
setUp()Set up all the internal data structures.
setWhichSingularTriplets(which)Set which singular triplets are to be sought.
solve()Solve the singular value problem.
valuesView([viewer])Display the computed singular values in a viewer.
vectorsView([viewer])Output computed singular vectors to a viewer.
view([viewer])Print the SVD data structure.
Attributes Summary
The direct solver (
DS) object associated.The maximum iteration count.
The type of the eigenvalue problem.
The tolerance.
Compute the residual norm of all approximate eigenpairs.
How to handle the transpose of the matrix.
The portion of the spectrum to be sought.
Methods Documentation
- appendOptionsPrefix(prefix=None)#
Append to the prefix used for searching for all SVD options in the database.
Logically collective.
- cancelMonitor()#
Clear all monitors for an
SVDobject.Logically collective.
See also
Source code at slepc4py/SLEPc/SVD.pyx:1149
- Return type:
- computeError(i, etype=None)#
Compute the error associated with the i-th singular triplet.
Collective.
Compute the error (based on the residual norm) associated with the i-th singular triplet.
- Parameters:
- Returns:
The error bound, computed in various ways from the residual norm \(\sqrt{\eta_1^2+\eta_2^2}\) where \(\eta_1 = \|A v - \sigma u\|_2\), \(\eta_2 = \|A^* u - \sigma v\|_2\), \(\sigma\) is the approximate singular value, \(u\) and \(v\) are the left and right singular vectors.
- Return type:
Notes
The index
ishould be a value between0andnconv-1(seegetConverged()).In the case of the GSVD, the two components of the residual norm are \(\eta_1 = \|s^2 A^*u-cB^*Bx\|_2\) and \(\eta_2 = ||c^2 B^*v-sA^*Ax||_2\), where \((\sigma,u,v,x)\) is the approximate generalized singular quadruple, with \(\sigma=c/s\).
See also
- create(comm=None)#
Create the SVD object.
Collective.
- Parameters:
comm (Comm | None) – MPI communicator; if not provided, it defaults to all processes.
- Return type:
See also
- destroy()#
Destroy the SVD object.
Collective.
See also
Source code at slepc4py/SLEPc/SVD.pyx:214
- Return type:
- errorView(etype=None, viewer=None)#
Display the errors associated with the computed solution.
Collective.
Display the errors and the singular values.
- 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 singular triplets and prints the singular values if all of them are below the requested tolerance. If the viewer has format
ASCII_INFO_DETAILthen a table with singular values and corresponding errors is printed.See also
- getBV()#
Get the basis vectors objects associated to the SVD object.
Not collective.
- Returns:
- Return type:
- getConverged()#
Get the number of converged singular triplets.
Not collective.
- Returns:
nconv – Number of converged singular triplets.
- Return type:
Notes
This function should be called after
solve()has finished.The value
nconvmay be different from the number of requested solutionsnsv, 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
- getCrossEPS()#
Get the eigensolver object associated to the singular value solver.
Collective.
- Returns:
The eigensolver object.
- Return type:
See also
- getCrossExplicitMatrix()#
Get the flag indicating if \(A^*A\) is built explicitly.
Not collective.
- Returns:
Trueif \(A^*A\) is built explicitly.- Return type:
- getCyclicEPS()#
Get the eigensolver object associated to the singular value solver.
Collective.
- Returns:
The eigensolver object.
- Return type:
See also
- getCyclicExplicitMatrix()#
Get the flag indicating if \(H(A)\) is built explicitly.
Not collective.
Get the flag indicating if \(H(A) = [ 0\; A ; A^T\; 0 ]\) is built explicitly.
- Returns:
Trueif \(H(A)\) is built explicitly.- Return type:
- getDS()#
Get the direct solver associated to the singular value solver.
Not collective.
- Returns:
The direct solver context.
- Return type:
- getDimensions()#
Get the number of singular values to compute and the dimension of the subspace.
Not collective.
- Returns:
- Return type:
See also
- getImplicitTranspose()#
Get the mode used to handle the transpose of the associated matrix.
Not collective.
- Returns:
How to handle the transpose (implicitly or not).
- 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
- getLanczosOneSide()#
Get if the variant of the Lanczos method to be used is one-sided or two-sided.
Not collective.
- Returns:
Trueif the method is one-sided.- Return type:
See also
- getMonitor()#
Get the list of monitor functions.
Not collective.
- Returns:
The list of monitor functions.
- Return type:
See also
- getOperators()#
Get the matrices associated with the singular value problem.
Collective.
- Returns:
A (
petsc4py.PETSc.Mat) – The matrix associated with the singular value problem.B (
petsc4py.PETSc.Mat) – The second matrix in the case of GSVD.
- Return type:
See also
- getOptionsPrefix()#
Get the prefix used for searching for all SVD options in the database.
Not collective.
- Returns:
The prefix string set for this SVD object.
- Return type:
- getProblemType()#
Get the problem type from the SVD object.
Not collective.
- Returns:
The problem type that was previously set.
- Return type:
See also
- getSignature(omega=None)#
Get the signature matrix defining a hyperbolic singular value problem.
Collective.
- Parameters:
omega (Vec | None) – Optional vector to store the diagonal elements of the signature matrix.
- Returns:
A vector containing the diagonal elements of the signature matrix.
- Return type:
See also
- getSingularTriplet(i, U=None, V=None)#
Get the i-th triplet of the singular value decomposition.
Collective.
Get the i-th triplet of the singular value decomposition as computed by
solve(). The solution consists of the singular value and its left and right singular vectors.- Parameters:
- Returns:
The computed singular value.
- Return type:
Notes
The index
ishould be a value between0andnconv-1(seegetConverged(). Singular triplets are indexed according to the ordering criterion established withsetWhichSingularTriplets().
- getStoppingTest()#
Get the stopping test function.
Not collective.
- Returns:
The stopping test function.
- Return type:
See also
- getTRLanczosExplicitMatrix()#
Get the flag indicating if \(Z=[A^*,B^*]^*\) is built explicitly.
Not collective.
- Returns:
Trueif \(Z=[A^*,B^*]^*\) is built explicitly.- Return type:
- getTRLanczosGBidiag()#
Get bidiagonalization choice used in the GSVD TRLanczos solver.
Not collective.
- Returns:
The bidiagonalization choice.
- Return type:
See also
- getTRLanczosKSP()#
Get the linear solver object associated with the SVD solver.
Collective.
- Returns:
The linear solver object.
- Return type:
See also
- getTRLanczosLocking()#
Get the locking flag used in the thick-restart Lanczos method.
Not collective.
- Returns:
The locking flag.
- Return type:
See also
- getTRLanczosOneSide()#
Get if the variant of the method to be used is one-sided or two-sided.
Not collective.
Get if the variant of the thick-restart Lanczos method to be used is one-sided or two-sided.
- Returns:
Trueif the method is one-sided.- Return type:
See also
- getTRLanczosRestart()#
Get the restart parameter used in the thick-restart Lanczos method.
Not collective.
- Returns:
The number of vectors to be kept at restart.
- Return type:
See also
- getThreshold()#
Get the threshold used in the threshold stopping test.
Not collective.
- Returns:
- Return type:
See also
- getTolerances()#
Get the tolerance and maximum iteration count.
Not collective.
Get the tolerance and maximum iteration count used by the default SVD convergence tests.
- Returns:
- Return type:
See also
- getTrackAll()#
Get the flag indicating if all residual norms must be computed or not.
Not collective.
- Returns:
Whether the solver computes all residuals or not.
- Return type:
See also
- getType()#
Get the SVD type of this object.
Not collective.
- Returns:
The solver currently being used.
- Return type:
See also
- getValue(i)#
Get the i-th singular value as computed by
solve().Collective.
- Parameters:
i (int) – Index of the solution to be obtained.
- Returns:
The computed singular value.
- Return type:
Notes
The index
ishould be a value between0andnconv-1(seegetConverged(). Singular triplets are indexed according to the ordering criterion established withsetWhichSingularTriplets().
- getVectors(i, U, V)#
Get the i-th left and right singular vectors as computed by
solve().Collective.
- Parameters:
- Return type:
Notes
The index
ishould be a value between0andnconv-1(seegetConverged(). Singular triplets are indexed according to the ordering criterion established withsetWhichSingularTriplets().
- getWhichSingularTriplets()#
Get which singular triplets are to be sought.
Not collective.
- Returns:
The singular values to be sought (either largest or smallest).
- Return type:
- isGeneralized()#
Tell if the SVD corresponds to a generalized singular value problem.
Not collective.
- Returns:
Trueif two matrices were set withsetOperators().- Return type:
See also
- isHyperbolic()#
Tell whether the SVD object corresponds to a hyperbolic singular value problem.
Not collective.
- Returns:
Trueif the problem was specified as hyperbolic.- Return type:
See also
- reset()#
Reset the SVD object.
Collective.
See also
Source code at slepc4py/SLEPc/SVD.pyx:228
- Return type:
- setBV(V, U=None)#
Set basis vectors objects associated to the SVD solver.
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:
See also
- setCrossEPS(eps)#
Set an eigensolver object associated to the singular value solver.
Collective.
See also
- setCrossExplicitMatrix(flag=True)#
Set if the eigensolver operator \(A^*A\) must be computed.
Logically collective.
Notes
In GSVD there are two cross product matrices, \(A^*A\) and \(B^*B\). In HSVD the expression for the cross product matrix is different, \(A^*\Omega A\).
By default the matrices are not built explicitly, but handled as shell matrices
- setCyclicEPS(eps)#
Set an eigensolver object associated to the singular value solver.
Collective.
See also
- setCyclicExplicitMatrix(flag=True)#
Set if the eigensolver operator \(H(A)\) must be computed explicitly.
Logically collective.
Set if the eigensolver operator \(H(A) = [ 0\; A ; A^T\; 0 ]\) must be computed explicitly.
Notes
In GSVD and HSVD the equivalent eigenvalue problem has generalized form, and hence two matrices are built.
By default the matrices are not built explicitly, but handled as shell matrices.
- setDS(ds)#
Set a direct solver object associated to the singular value solver.
Collective.
- setDimensions(nsv=None, ncv=None, mpd=None)#
Set the number of singular values to compute and the dimension of the subspace.
Logically collective.
- 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
nsvis small, the user setsncv(a reasonable default is 2 *nsv).In cases where
nsvis large, the user setsmpd.
The value of
ncvshould always be betweennsvand (nsv+mpd), typicallyncv=nsv+mpd. Ifnsvis not too large,mpd=nsvis a reasonable choice, otherwise a smaller value should be used.See also
- setFromOptions()#
Set SVD 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/SVD.pyx:376
- Return type:
- setImplicitTranspose(mode)#
Set how to handle the transpose of the associated matrix.
Logically collective.
Notes
By default, the transpose of the matrix is explicitly built (if the matrix has defined the
Mat.transpose()operation).If this flag is set to
True, the solver does not build the transpose, but handles it implicitly viaMat.multTranspose()(orMat.multHermitianTranspose()in the complex case).See also
- setInitialSpace(spaceright=None, spaceleft=None)#
Set the initial spaces from which the SVD solver starts to iterate.
Collective.
- Parameters:
- Return type:
Notes
The initial right and left spaces are rough approximations to the right and/or left singular subspaces from which the solver starts to iterate. It is not necessary to provide both sets of vectors.
Some solvers start to iterate on a single vector (initial vector). In that case, the other vectors are ignored.
These vectors do not persist from one
solve()call to the other, so the initial spaces 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 singular spaces. Then, convergence may be faster.
See also
- setLanczosOneSide(flag=True)#
Set if the variant of the Lanczos method to be used is one-sided or two-sided.
Logically collective.
Notes
By default, a two-sided variant is selected, which is sometimes slightly more robust. However, the one-sided variant is faster because it avoids the orthogonalization associated to left singular vectors. It also saves the memory required for storing such vectors.
See also
- setMonitor(monitor, args=None, kargs=None)#
Append a monitor function to the list of monitors.
Logically collective.
See also
- setOperators(A, B=None)#
Set the matrices associated with the singular value problem.
Collective.
- Parameters:
- Return type:
See also
- setOptionsPrefix(prefix=None)#
Set the prefix used for searching for all SVD options in the database.
Logically collective.
- Parameters:
prefix (str | None) – The prefix string to prepend to all SVD 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 SVD contexts, one could call:
S1.setOptionsPrefix("svd1_") S2.setOptionsPrefix("svd2_")
- setProblemType(problem_type)#
Set the type of the singular value problem.
Logically collective.
- Parameters:
problem_type (ProblemType) – The problem type to be set.
- Return type:
Notes
The GSVD requires that two matrices have been passed via
setOperators(). The HSVD requires that a signature matrix has been passed viasetSignature().See also
setOperators,setSignature,getProblemType,SVDSetProblemType
- setSignature(omega=None)#
Set the signature matrix defining a hyperbolic singular value problem.
Collective.
- Parameters:
omega (Vec | None) – A vector containing the diagonal elements of the signature matrix.
- Return type:
See also
- setStoppingTest(stopping, args=None, kargs=None)#
Set a function to decide when to stop the outer iteration of the eigensolver.
Logically collective.
See also
- setTRLanczosExplicitMatrix(flag=True)#
Set if the matrix \(Z=[A^*,B^*]^*\) must be built explicitly.
Logically collective.
Notes
This option is relevant for the GSVD case only. \(Z\) is the coefficient matrix of the least-squares solver used internally.
- setTRLanczosGBidiag(bidiag)#
Set the bidiagonalization choice to use in the GSVD TRLanczos solver.
Logically collective.
- Parameters:
bidiag (TRLanczosGBidiag) – The bidiagonalization choice.
- Return type:
See also
- setTRLanczosKSP(ksp)#
Set a linear solver object associated to the SVD solver.
Collective.
See also
- setTRLanczosLocking(lock)#
Toggle between locking and non-locking variants of TRLanczos.
Logically collective.
Notes
The default is to lock converged singular triplets 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
- setTRLanczosOneSide(flag=True)#
Set if the variant of the method to be used is one-sided or two-sided.
Logically collective.
Set if the variant of the thick-restart Lanczos method to be used is one-sided or two-sided.
Notes
By default, a two-sided variant is selected, which is sometimes slightly more robust. However, the one-sided variant is faster because it avoids the orthogonalization associated to left singular vectors.
See also
- setTRLanczosRestart(keep)#
Set the restart parameter for the thick-restart Lanczos method.
Logically collective.
Set the restart parameter for the thick-restart Lanczos method, in particular 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
- setThreshold(thres, rel=False)#
Set the threshold used in the threshold stopping test.
Logically collective.
- Parameters:
- Return type:
Notes
This function internally sets a special stopping test based on the threshold, where singular values are computed in sequence until one of the computed singular values is below/above the threshold (depending on whether largest or smallest singular values are computed).
In the case of largest singular values, the threshold can be made relative with respect to the largest singular value (i.e., the matrix norm).
The details are given in
SVDSetThreshold.See also
- setTolerances(tol=None, max_it=None)#
Set the tolerance and maximum iteration count used.
Logically collective.
Set the tolerance and maximum iteration count used by the default SVD convergence tests.
- Parameters:
- Return type:
Notes
Use
DETERMINEformax_itto assign a reasonably good value, which is dependent on the solution method.See also
- setTrackAll(trackall)#
Set flag to compute the residual of all singular triplets.
Logically collective.
Set if the solver must compute the residual of all approximate singular triplets or not.
See also
- setType(svd_type)#
Set the particular solver to be used in the SVD object.
Logically collective.
Notes
The default is
CROSS. Normally, it is best to usesetFromOptions()and then set the SVD 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 singular value solver.
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/SVD.pyx:1164
- Return type:
- setWhichSingularTriplets(which)#
Set which singular triplets are to be sought.
Logically collective.
- Parameters:
which (Which) – The singular values to be sought (either largest or smallest).
- Return type:
- solve()#
Solve the singular value problem.
Collective.
Notes
The problem matrices are specified with
setOperators().solve()will return without generating an error regardless of whether all requested solutions were computed or not. CallgetConverged()to get the actual number of computed solutions, andgetConvergedReason()to determine if the solver converged or failed and why.See also
setUp,setOperators,getConverged,getConvergedReason,SVDSolveSource code at slepc4py/SLEPc/SVD.pyx:1185
- Return type:
- valuesView(viewer=None)#
Display the computed singular values 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 singular vectors 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 SVD 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.
- problem_type#
The type of the eigenvalue problem.
- tol#
The tolerance.
- track_all#
Compute the residual norm of all approximate eigenpairs.
- transpose_mode#
How to handle the transpose of the matrix.
- which#
The portion of the spectrum to be sought.