NAME
PDL::SVDLIBC - PDL interface to Doug Rohde's SVD C Library
SYNOPSIS
use PDL;
use PDL::SVDLIBC;
##---------------------------------------------------------------------
## Input matrix (dense)
##---------------------------------------------------------------------
$n = 100; ##-- number of columns
$m = 50; ##-- number of rows
$a = random(double,$n,$m); ##-- random matrix
##---------------------------------------------------------------------
## Output pdls
##---------------------------------------------------------------------
$d = $n; ##-- max number of output dimensions
$ut = zeroes(double,$m,$d); ##-- left singular components
$s = zeroes(double,$d); ##-- singular values (diagnonal vector)
$vt = zeroes(double,$n,$n); ##-- right singular components
##---------------------------------------------------------------------
## Singular Value Decomposition (dense)
##---------------------------------------------------------------------
svdlas2d($a, $maxiters, $end, $kappa, $ut, $s, $vt);
##---------------------------------------------------------------------
## Singular Value Decomposition (sparse)
##---------------------------------------------------------------------
use PDL::CCS;
($ptr,$rowids,$nzvals) = ccsencode($a);
$ptr->reshape($ptr->nelem+1);
$ptr->set(-1, $rowids->nelem);
svdlas2($ptr, $rowids, $nzvals, $m, $maxiters, $end, $kappa, $ut, $s, $vt);
DESCRIPTION
PDL::SVDLIBC provides a PDL interface to the SVDLIBC routines for singular value decomposition of large sparse matrices. SVDLIBC is available from http://tedlab.mit.edu/~dr/SVDLIBC/
FUNCTIONS
SVDLIBC Globals
There are several global data structures still lurking in the SVDLIBC code, so expect problems if you are trying to run more than one 'las2' procedure at once (even in different processes).
PDL::SVDLIBC provides access to (some of) the SVDLIBC globals through the following functions, which are not exported.
PDL::SVDLIBC::verbosity()
PDL::SVDLIBC::verbosity($level)
Get/set the current SVDLIBC verbosity level. Valid values for $level are between 0 (no messages) and 2 (many messages).
PDL::SVDLIBC::svdVersion()
Returns a string representing the SVDLIBC version this module was compiled with.
SVD Utilities
_svdccsencode
Signature: (double a(n,m); indx [o]ptr(n1); indx [o]rowids(nnz); double [o]nzvals(nnz))
info not available
_svdccsencode does not process bad values. It will set the bad-value flag of all output piddles if the flag is set for any of the input piddles.
svdlas2a
indx ptr(nplus1);
indx rowids(nnz);
double nzvals(nnz);
indx nrows(); ##-- default: max($rowids)+1
int d(); ##-- default: nplus1-1
int iterations(); ##-- default: 2*$d
double end(2); ##-- default: [-1e-30,1e-30]
double kappa(); ##-- default: 1e-6
double [o]ut(m,d); ##-- default: new
double [o] s(d); ##-- default: new
double [o]vt(n,d); ##-- default: new
Uses a variant of the single-vector Lanczos method (Lanczos, 1950) to compute the singular value decomposition of a sparse matrix with $nrows() rows and data encoded in Harwell-Boeing sparse format in the input parameters $ptr(), $rowids(), and $nzvals(). See "PDL::CCS" for a way to acquire these parameters from a dense input matrix, but note that for svdlas2(), the column pointer $ptr() is of size ($n+1) for a dense matrix $a with $n columns, where $ptr($n)==$nnz is the total number of nonzero values in $a.
$iterations() is the maximum number of Lanczos iterations to perform.
$end() specifies two endpoints of an interval within which all unwanted eigenvalues lie.
$kappa() is a double containing the relative accuracy of Ritz values acceptable as eigenvalues.
The left singular components are returned in the matrix $ut(), the singular values themselved in the vector $s(), and the right singular components in the matrix $vt(). Note that $ut() and $vt() are transposed, and must be specified explicitly in the call, so that the degree of reduction (the size parameter $d) can be determined. If $d==$n, then a full decomposition will be computed, and on return, $ut() and $vt() should be transposed instances of the matrices $u() and $v() as returned by PDL::MatrixOps::svd().
The Lanczos method as used here seems to be consistently the fastest. This algorithm has the drawback that the low order singular values may be relatively imprecise, but that is not a problem for most users who only want the higher-order values or who can tolerate some imprecision.
See also: svdlas2d()
svdlas2
Signature: (
indx ptr(nplus1);
indx rowids(nnz);
double nzvals(nnz);
indx nrows();
int iterations();
double end(2);
double kappa();
double [o]ut(m,d);
double [o] s(d);
double [o]vt(n,d);
)
Guts for svdlas2a(). No default instantiation, and slightly different calling conventions.
svdlas2 does not process bad values. It will set the bad-value flag of all output piddles if the flag is set for any of the input piddles.
svdlas2ad
double a(n,m);
int d(); ##-- default: $n
int iterations(); ##-- default: 2*$d
double end(2); ##-- default: [-1e-30,1e-30]
double kappa(); ##-- default: 1e-6
double [o]ut(m,d); ##-- default: new
double [o] s(d); ##-- default: new
double [o]vt(n,d); ##-- default: new
As for svdlas2(), but implicitly converts the dense input matrix $a() to sparse format before computing the decomposition.
svdlas2d
Signature: (
double a(n,m);
int iterations();
double end(2);
double kappa();
double [o]ut(m,d);
double [o] s(d);
double [o]vt(n,d);
)
Guts for _svdlas2d().
svdlas2d does not process bad values. It will set the bad-value flag of all output piddles if the flag is set for any of the input piddles.
ACKNOWLEDGEMENTS
Perl by Larry Wall.
PDL by Karl Glazebrook, Tuomas J. Lukka, Christian Soeller, and others.
SVDLIBC by Dough Rohde.
SVDPACKC by Michael Berry, Theresa Do, Gavin O'Brien, Vijay Krishna and Sowmini Varadhan.
KNOWN BUGS
Globals still lurk in the depths of SVDLIBC.
AUTHOR
Bryan Jurish <moocow@cpan.org>
Copyright Policy
Copyright (C) 2005-2013, Bryan Jurish. All rights reserved.
This package is free software, and entirely without warranty. You may redistribute it and/or modify it under the same terms as Perl itself.
SEE ALSO
perl(1), PDL(3perl), PDL::CCS(3perl), SVDLIBC documentation.