slepc-main 2024-12-17
Report Typos and Errors

SVDGetSingularTriplet

Gets the i-th triplet of the singular value decomposition as computed by SVDSolve(). The solution consists in the singular value and its left and right singular vectors.

Synopsis

#include "slepcsvd.h" 
PetscErrorCode SVDGetSingularTriplet(SVD svd,PetscInt i,PetscReal *sigma,Vec u,Vec v)
Collective

Input Parameters

svd  - singular value solver context
i  - index of the solution

Output Parameters

sigma  - singular value
u  - left singular vector
v  - right singular vector

Note

Both u or v can be NULL if singular vectors are not required. Otherwise, the caller must provide valid Vec objects, i.e., they must be created by the calling program with e.g. MatCreateVecs().

The index i should be a value between 0 and nconv-1 (see SVDGetConverged()). Singular triplets are indexed according to the ordering criterion established with SVDSetWhichSingularTriplets().

In the case of GSVD, the solution consists in three vectors u,v,x that are returned as follows. Vector x is returned in the right singular vector (argument v) and has length equal to the number of columns of A and B. The other two vectors are returned stacked on top of each other [u;v] in the left singular vector argument, with length equal to m+n (number of rows of A plus number of rows of B).

See Also

SVDSolve(), SVDGetConverged(), SVDSetWhichSingularTriplets()

Level

beginner

Location

src/svd/interface/svdsolve.c

Examples

src/svd/tutorials/ex8.c
src/svd/tutorials/ex15.c
src/svd/tutorials/ex15f.F90
src/svd/tutorials/ex45.c
src/svd/tutorials/ex52.c
src/svd/tutorials/ex53.c


Index of all SVD routines
Table of Contents for all manual pages
Index of all manual pages