Rcpp_matrix_vector_blockMult

C++ Function Reference

1 Signature

Rcpp::RObject BigDataStatMeth::Rcpp_matrix_vector_blockMult(T A, T B, Rcpp::Nullable< bool > bparal, Rcpp::Nullable< int > iblock_size, Rcpp::Nullable< int > threads)

2 Description

Block-based matrix-vector multiplication.

3 Parameters

  • A (T): Input matrix
  • B (T): Input vector
  • bparal (Rcpp::Nullable< bool >): Whether to use parallel processing
  • iblock_size (Rcpp::Nullable< int >): Block size for computation
  • threads (Rcpp::Nullable< int >): Number of threads for parallel processing

4 Returns

Result of matrix-vector multiplication

5 Details

Implements block-based matrix-vector multiplication with:Optional parallel processingConfigurable block sizesThread count controlCache-friendly algorithms

6 Call Graph

Function dependencies

7 Source Code

File: inst/include/memAlgebra/memMultiplication.hppLines 394-494

inline Rcpp::RObject Rcpp_matrix_vector_blockMult( T  A, T  B, Rcpp::Nullable<bool> bparal, 
                            Rcpp::Nullable<int> iblock_size, Rcpp::Nullable<int> threads)
    {
        
        // NOTA: Per defecte, multiplica per columnes tal i com raja.... 

        bool btransposed = false;
        // unsigned int ithreads;
        hsize_t block_size;
        // int chunks;
        
        Rcpp::NumericMatrix X = Rcpp::as<Rcpp::NumericMatrix>(A);
        Rcpp::NumericVector Y = Rcpp::as<Rcpp::NumericVector>(B);
        Rcpp::NumericMatrix C;
        
        // Matrix
        hsize_t M = X.rows(), N = X.cols();
        // Vector
        hsize_t K = Y.length();
        
        try {
            
            if( K==N || K==M) {
                if ( K == N){
                    // multiplies vector to every col
                    btransposed = true;
                    
                    X = Rcpp::transpose(X);
                    
                    N = X.rows();
                    M = X.cols();
                } 
                
                std::vector<hsize_t> vsizetoRead;
                std::vector<hsize_t> vstart;
                
                // ithreads = get_number_threads(threads, bparal);
                
                C = Rcpp::no_init( M, N);
                
                if( iblock_size.isNotNull()) {
                    block_size =  Rcpp::as<int>(iblock_size);  
                } else {
                    block_size = getMatrixBlockSize( N, M).at(0);
                }
                
                // minimum block size: 2 columns
                if(block_size <= 0 ) {
                    block_size = M*2;
                }
                
                // Mínimum block size: 2 columns
                getBlockPositionsSizes( M*N, block_size, vstart, vsizetoRead );
                
                // chunks = vstart.size()/ithreads;
                
                #pragma omp parallel num_threads( get_number_threads(threads, bparal) ) shared(A, B, C) //, chunks)
                {
                #pragma omp for schedule (static)
                    for (hsize_t ii = 0; ii < vstart.size(); ii ++)
                    {
                        // Duplicate vector
                        std::size_t const no_of_duplicates = vsizetoRead[ii] / Y.length();
                        
                        std::vector<double> v = Rcpp::as<std::vector<double> >(Y); 
                        v.reserve(Y.size() * no_of_duplicates);
                        auto end = std::end(v);
                        
                        for(std::size_t i = 1; i < no_of_duplicates; ++i)
                            v.insert(std::end(v), std::begin(v), end);
                        
                        // Mult vector to matrix by columns / rows
                        if( vstart[ii] + vsizetoRead[ii] >= M*N ) {
                            std::transform (X.begin() + vstart[ii], X.end(),
                                            v.begin(), C.begin() + vstart[ii], std::multiplies<double>());
                        } else {
                            std::transform (X.begin() + vstart[ii], X.begin() + vstart[ii] + vsizetoRead[ii],
                                            v.begin() , C.begin() + vstart[ii], std::multiplies<double>());   
                        }
                    }
                }

            } else {
                
                Rcpp::Rcout<< "vector sum error: non-conformable arguments\n";
                return(R_NilValue);
            }
            
        } catch(std::exception& ex) {
            Rcpp::Rcout<< "c++ exception Rcpp_matrix_vector_blockMult: "<<ex.what()<< " \n";
            return(R_NilValue);
        }
        
        if(btransposed == true){
            Rcpp::transpose(C);
        } 
        
        C.attr("dim") = Rcpp::Dimension( M, N);
        
        return(C);
    }

8 Usage Example

#include "BigDataStatMeth.hpp"

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