The Perl Toolchain Summit needs more sponsors. If your company depends on Perl, please support this very important event.

NAME

PDL::SVDSLEPc - PDL interface to SLEPc sparse singular value decomposition

SYNOPSIS

 use PDL;
 use PDL::SVDSLEPc;

 ##---------------------------------------------------------------------
 ## Input matrix (dense)
 
 $m = 100;                    ##-- number of rows
 $n = 50;                     ##-- number of columns
 $a = random(double,$m,$n);   ##-- random matrix

 ##---------------------------------------------------------------------
 ## Input matrix (sparse)
 
 use PDL::CCS;
 $a  *= ($a->random>0.9);      ##-- make sparse
 $ccs = $a->toccs();           ##-- encode as PDL::CCS::Nd object
 
 ##-- get Harwell-Boeing encoding triple
 $ptr    = $ccs->ptr(0);
 $colids = $ccs->_whichND->slice("(1),");
 $nzvals = $cccs->_nzvals;

 ##---------------------------------------------------------------------
 ## SLEPc Singular Value Decomposition
 
 slepc_svd_help();            ##-- print available options to STDOUT
 
 ($u,$s,$v) = slepc_svd($ccs, ['-svd_nsv'=>32]);                 ##-- from PDL::CCS::Nd object
 ($u,$s,$v) = slepc_svd($ptr,$colids,$nzvals, ['-svd_nsv'=>32]); ##-- from Harwell-Boeing encoding

DESCRIPTION

PDL::SVDSLEPc provides a PDL interface to the SLEPc singular value decomposition solver(s). SLEPc itself is available from http://slepc.upv.es/.

FUNCTIONS

CONSTANTS

PDL::SVDSLEPc provides access to the following SLEPc constants:

PDL::SVDSLEPc::slepc_version()

Returns a string representing the SLEPc library version this module was compiled with.

PDL::SVDSLEPc::petsc_version()

Returns a string representing the PETSc library version this module was compiled with.

PDL::SVDSLEPc::library_version()

In list context returns a pair (slepc_version(),petsc_version()), in scalar context returns a string with both versions separated with a semicolon.

PDL::SVDSLEPc::MPI_Comm_size()

Returns the number of MPI processes available. Using multiple MPI processes with mpiexec behaves strangely with perl at the moment, so this should generally return 1.

SVD UTILITIES

The following functions are provided as quick and dirty wrappers for the SLEPc SVD solver class.

slepc_svd_help

Prints a help message with all supported SLEPc SVD options to STDOUT. Really just a wrapper for

 slepc_svd(null,null,null, 0,0,0, ['-help']);

slepc_svd

  Signature: (
              int  rowptr(mplus1);
              int  colids(nnz);
              double nzvals(nnz);
              double [o]u(d,n);         # optional
              double [o]s(d);           # optional
              double [o]v(d,m);         # optional
              int    M=>d;              # optional
              int    N=>n;              # optional
              int    D=>d;              # optional
              \@options;                # optional
             )

Compute the (truncated) singular value decomposition of a sparse matrix using a SLEPc SVD solver. The sparse input matrix with logical (row-primary) dimensions ($m,$n) is passed in encoded in Harwell-Boeing sparse format in the input parameters $rowptr(), $colids(), and $nzvals(). See PDL::CCS for a way to acquire these parameters from a dense input matrix, but note that for this function, the row-pointer $rowptr() is of size ($m+1) for a dense matrix $a with $m rows, where $rowptr($m)==$nnz is the total number of nonzero values in $a. As an alternative, a single PDL::CCS::Nd object can be passed in place of of $rowptr, $colids, and $nzvals.

The left singular components are returned in the matrix $u(), the singular values themselved in the vector $s(), and the right singular components in the matrix $v(). These output piddles, as well as the logical dimensions ($m,$n) of the input matrix and the size $d of the truncated SVD to be computed may be specified explicitly in the call, but otherwise will be estimated from $rowptr(), $colids(), and/or \@options. If $d==min($m,$n) [the default], then a full decomposition will be computed, and on return, $u() and $v() should be variants (up to sign and specified error tolerance) of the matrices returned by PDL::MatrixOps::svd().

Additional options to the underlying SLEPc and PETSc code can be passed command-line style in the ARRAY-ref \@options; see the output of slepc_svd_help() for a list of available options. In particular, the option -svd_nsv can be used to specify the number of singular values to be returned ($d) if you choose to omit both the $d paramter and nontrivial output piddles. For example, in order to compute a truncated SVD using with 32 singular values using the Thick-restart Lanczos method with at most 128 iterations and a tolerance of 1e-5, you could call:

 ($u,$s,$v) = slepc_svd($rowptr,$colids,$nzvals,
                        [qw(-svd_type trlanczos
                            -svd_nsv 32
                            -svd_max_it 128
                            -svd_tol 1e-5)]);

... or if you already have a PDL::CCS::Nd object $ccs handy:

  ($u,$s,$v) = $ccs->slepc_svd([-svd_type=>'trlanczos', -svd_nsv=>32, -svd_max_it=>128, -svd_tol=>1e-5]);

_slepc_svd_int

  Signature: (
    int  rowptr(mplus1);
    int  colids(nnz);
    double nzvals(nnz);
    double [o]u(d,n);
    double [o]s(d);
    double [o]v(d,m);
    ; 
    int M=>m;
    int N=>n;
    int D=>d;
    IV optsArray;
    )

Guts for slepc_svd() with stricter calling conventions: The input matrix must be passed as a Harwell-Boeing triple ($rowptr,$colids,$nzvals), and the size parameters M, N, and D and options array optsArray are all mandatory.

_slepc_svd_int 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.

SLEPc by Carmen Campos, Jose E. Roman, Eloy Romero, and Andres Tomas.

KNOWN BUGS

no abstract interface

There should really be a more general and abstract PDL interface to SLEPc/PETsc.

OpenMPI Errors "mca: base: component find: unable to open ..."

You might see OpenMPI errors such as the following when trying to use this module:

 mca: base: component find: unable to open /usr/lib/openmpi/lib/openmpi/mca_paffinity_hwloc: perhaps a missing symbol, or compiled for a different version of Open MPI? (ignored)

If you do, you probably need to configure your runtime linker to pre-load the OpenMPI libraries, e.g. with

 export LD_PRELOAD=/usr/lib/libmpi.so

or similar. An alternative is to build OpenMPI with the --disable-dlopen option. See http://www.open-mpi.org/faq/?category=troubleshooting#missing-symbols for details.

OpenMPI warnings "... unable to find any relevant network interfaces ... (openib)"

This OpenMPI warning has been observed on Ubuntu 14.04; it can be suppressed by setting the OpenMPI MCA btl ("byte transfer layer") parameter to exclude the openib module. This can be accomplished in various ways, e.g.:

via command-line parameters to mpiexec:

Call your program as:

 $ mpiexec --mca btl ^openib PROGRAM...
via environment variables

You can set the OpenMPI MCA paramters via environment variables, e.g.:

 $ export OMPI_MCA_btl="^openib"
 $ PROGRAM...
via configuration files

You can set OpenMPI MCA parameters via $HOME/.openmpi/mac-params.conf:

 ##-- suppress annoying warnings about missing openib
 btl = ^openib

See http://www.open-mpi.de/faq/?category=tuning#setting-mca-params for more details.

AUTHOR

Bryan Jurish <moocow@cpan.org>

COPYRIGHT AND LICENSE

Copyright (c) 2015, 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, either version 5.20.2 or any newer version of Perl 5 you have available.

SEE ALSO

perl(1), PDL(3perl), PDL::CCS(3perl), PDL::SVDLIBC(3perl), the SLEPc documentation at http://slepc.upv.es/documentation/current/docs/index.html.