mlpack
Public Member Functions | List of all members
mlpack::svd::RandomizedBlockKrylovSVD Class Reference

Randomized block krylov SVD is a matrix factorization that is based on randomized matrix approximation techniques, developed in in "Randomized Block Krylov Methods for Stronger and Faster Approximate Singular Value Decomposition". More...

#include <randomized_block_krylov_svd.hpp>

Public Member Functions

 RandomizedBlockKrylovSVD (const arma::mat &data, arma::mat &u, arma::vec &s, arma::mat &v, const size_t maxIterations=2, const size_t rank=0, const size_t blockSize=0)
 Create object for the randomized block krylov SVD method. More...
 
 RandomizedBlockKrylovSVD (const size_t maxIterations=2, const size_t blockSize=0)
 Create object for the randomized block krylov SVD method. More...
 
void Apply (const arma::mat &data, arma::mat &u, arma::vec &s, arma::mat &v, const size_t rank)
 Apply Principal Component Analysis to the provided data set using the randomized block krylov SVD. More...
 
size_t MaxIterations () const
 Get the number of iterations for the power method.
 
size_t & MaxIterations ()
 Modify the number of iterations for the power method.
 
size_t BlockSize () const
 Get the block size.
 
size_t & BlockSize ()
 Modify the block size.
 

Detailed Description

Randomized block krylov SVD is a matrix factorization that is based on randomized matrix approximation techniques, developed in in "Randomized Block Krylov Methods for Stronger and Faster Approximate Singular Value Decomposition".

For more information, see the following.

@inproceedings{Musco2015,
author = {Cameron Musco and Christopher Musco},
title = {Randomized Block Krylov Methods for Stronger and Faster
Approximate Singular Value Decomposition},
booktitle = {Advances in Neural Information Processing Systems 28: Annual
Conference on Neural Information Processing Systems 2015,
December 7-12, 2015, Montreal, Quebec, Canada},
pages = {1396--1404},
year = {2015},
}

An example of how to use the interface is shown below:

arma::mat data; // Rating data in the form of coordinate list.
const size_t rank = 20; // Rank used for the decomposition.
// Make a RandomizedBlockKrylovSVD object.
arma::mat u, s, v;
// Use the Apply() method to get a factorization.
bSVD.Apply(data, u, s, v, rank);

Constructor & Destructor Documentation

◆ RandomizedBlockKrylovSVD() [1/2]

mlpack::svd::RandomizedBlockKrylovSVD::RandomizedBlockKrylovSVD ( const arma::mat &  data,
arma::mat &  u,
arma::vec &  s,
arma::mat &  v,
const size_t  maxIterations = 2,
const size_t  rank = 0,
const size_t  blockSize = 0 
)

Create object for the randomized block krylov SVD method.

Parameters
dataData matrix.
uFirst unitary matrix.
vSecond unitary matrix.
sDiagonal matrix of singular values.
maxIterationsNumber of iterations for the power method (Default: 2).
rankRank of the approximation (Default: number of rows.)
blockSizeThe block size, must be >= rank (Default: rank + 10).

◆ RandomizedBlockKrylovSVD() [2/2]

mlpack::svd::RandomizedBlockKrylovSVD::RandomizedBlockKrylovSVD ( const size_t  maxIterations = 2,
const size_t  blockSize = 0 
)

Create object for the randomized block krylov SVD method.

Parameters
maxIterationsNumber of iterations for the power method (Default: 2).
blockSizeThe block size, must be >= rank (Default: rank + 10).

Member Function Documentation

◆ Apply()

void mlpack::svd::RandomizedBlockKrylovSVD::Apply ( const arma::mat &  data,
arma::mat &  u,
arma::vec &  s,
arma::mat &  v,
const size_t  rank 
)

Apply Principal Component Analysis to the provided data set using the randomized block krylov SVD.

Parameters
dataData matrix.
uFirst unitary matrix.
vSecond unitary matrix.
sDiagonal matrix of singular values.
rankRank of the approximation.

The documentation for this class was generated from the following files: