RcppbdSVD

C++ Function Reference

1 Signature

svdeig BigDataStatMeth::RcppbdSVD(Eigen::MatrixXd &X, int k, int ncv, bool bcenter, bool bscale)

2 Description

Compute SVD decomposition using Spectra eigenvalue solver.

3 Parameters

  • X (Eigen::MatrixXd &): Input matrix for SVD computation
  • k (int): Number of singular values/vectors to compute (0 = auto-select)
  • ncv (int): Number of Arnoldi vectors to use (0 = auto-select)
  • bcenter (bool): Whether to center the data before computation
  • bscale (bool): Whether to scale the data before computation

4 Returns

svdeig Structure containing U, S, V matrices and computation status

5 Details

Performs Singular Value Decomposition of a matrix using the Spectra library for eigenvalue computation. This function computes SVD by solving the eigenvalue problems XTX and XXT.

6 Call Graph

Function dependencies

7 Source Code

File: inst/include/hdf5Algebra/matrixSvd.hppLines 68-136

inline svdeig RcppbdSVD( Eigen::MatrixXd& X, int k, int ncv, bool bcenter, bool bscale )
    {
        
        svdeig retsvd;
        Eigen::MatrixXd nX;
        int nconv [[maybe_unused]];
        
        if( k==0 )    k = (std::min(X.rows(), X.cols()))-1;
        else if (k > (std::min(X.rows(), X.cols()))-1 ) k = (std::min(X.rows(), X.cols()))-1;
        
        if(ncv == 0)  ncv = k + 1 ;
        if(ncv<k) ncv = k + 1;
        
        {
            Eigen::MatrixXd Xtcp;
            if(bcenter ==true || bscale == true)  {
                nX = RcppNormalize_Data(X, bcenter, bscale, false);
                Xtcp =  bdtcrossproduct(nX);
            } else {
                Xtcp =  bdtcrossproduct(X);
            }
            
            Spectra::DenseSymMatProd<double> op(Xtcp);
            // Updated for Spectra 1.0.1: removed template parameters, pass object by reference
            Spectra::SymEigsSolver<Spectra::DenseSymMatProd<double>> eigs(op, k, ncv);
            
            // Initialize and compute
            eigs.init();
            // Updated for Spectra 1.0.1: SortRule as runtime parameter
            nconv = eigs.compute(Spectra::SortRule::LargestAlge);
            
            // Updated for Spectra 1.0.1: enum class for status check
            if(eigs.info() == Spectra::CompInfo::Successful) {
                retsvd.d = eigs.eigenvalues().cwiseSqrt();
                retsvd.u = eigs.eigenvectors();
                retsvd.bokuv = true;
            } else {
                retsvd.bokuv = false;
            }
        }
        if(retsvd.bokuv == true)
        {
            Eigen::MatrixXd Xcp;
            if(bcenter ==true || bscale==true ) {
                Xcp =  bdcrossproduct(nX);  
            } else {
                Xcp =  bdcrossproduct(X);
            }  
            
            Spectra::DenseSymMatProd<double> opv(Xcp);
            // Updated for Spectra 1.0.1: removed template parameters, pass object by reference
            Spectra::SymEigsSolver<Spectra::DenseSymMatProd<double>> eigsv(opv, k, ncv);
            
            // Initialize and compute
            eigsv.init();
            // Updated for Spectra 1.0.1: SortRule as runtime parameter
            nconv = eigsv.compute(Spectra::SortRule::LargestAlge);
            
            // Retrieve results
            // Updated for Spectra 1.0.1: enum class for status check
            if(eigsv.info() == Spectra::CompInfo::Successful) {
                retsvd.v = eigsv.eigenvectors();
            } else {
                retsvd.bokd = false;
            }
        }
        
        return retsvd;
    }

8 Usage Example

#include "BigDataStatMeth.hpp"

// Example usage
auto result = RcppbdSVD(...);