slepc-3.13.3 2020-06-14
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)
Not Collective, but vectors are shared by all processors that share the SVD
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().
See Also
SVDSolve(), SVDGetConverged(), SVDSetWhichSingularTriplets()
Location: src/svd/interface/svdsolve.c
Examples
src/svd/tutorials/ex8.c.html
src/svd/tutorials/ex15.c.html
src/svd/tutorials/ex15f.F.html
Index of all SVD routines
Table of Contents for all manual pages
Index of all manual pages