BVSVDAndRank#
Compute the SVD (left singular vectors only, and singular values) and determine the numerical rank according to a tolerance.
Synopsis#
#include "slepcbv.h"
PetscErrorCode BVSVDAndRank(BV S,PetscInt m,PetscInt l,PetscReal delta,BVSVDMethod meth,PetscScalar *A,PetscReal *sigma,PetscInt *rank)
Collective
Input Parameters#
S - the basis vectors
m - the moment degree
l - the block size
delta - the tolerance used to determine the rank
meth - the method to be used
Output Parameters#
A - workspace, on output contains relevant values in the CAA method
sigma - computed singular values
rank - estimated rank (optional)
Notes#
This function computes [U,Sigma,V] = svd(S) and replaces S with U. The current implementation computes this via S’*S, and it may include some kind of iterative refinement to improve accuracy in some cases.
The parameters m and l refer to the moment and block size of contour integral methods. All columns up to ml are modified, and the active columns are set to 0..ml.
The method is one of BV_SVD_METHOD_REFINE, BV_SVD_METHOD_QR, BV_SVD_METHOD_QR_CAA.
The A workspace should be mlm*l in size.
Once the decomposition is computed, the numerical rank is estimated by counting the number of singular values that are larger than the tolerance delta, relative to the first singular value.
See Also#
BVSetActiveColumns()
Level#
developer
Location#
src/sys/classes/bv/interface/bvcontour.c
Index of all BV routines Table of Contents for all manual pages Index of all manual pages