Bblock_matrix_mul_parallel
C++ Function Reference
1 Signature
Eigen::MatrixXd BigDataStatMeth::Bblock_matrix_mul_parallel(Eigen::MatrixXd A, Eigen::MatrixXd B, int block_size, Rcpp::Nullable< int > threads=R_NilValue)2 Description
Parallel block-based matrix multiplication.
3 Parameters
A(Eigen::MatrixXd): First input matrixB(Eigen::MatrixXd): Second input matrixblock_size(int): Size of blocks for computationthreads(Rcpp::Nullable< int >): Number of threads for parallel processing
4 Returns
Result of matrix multiplication
5 Details
Implements parallel block-based matrix multiplication for in-memory matrices.
6 Call Graph
7 Source Code
NoteImplementation
File: inst/include/hdf5Algebra/multiplication 2.hpp • Lines 110-171
inline Eigen::MatrixXd Bblock_matrix_mul_parallel( Eigen::MatrixXd A, Eigen::MatrixXd B,
int block_size, Rcpp::Nullable<int> threads = R_NilValue)
{
// unsigned int ithreads;
Eigen::MatrixXd C;
try {
int M = A.rows();
int K = A.cols();
int N = B.cols();
C = Eigen::MatrixXd::Zero(M,N) ;
if( A.rows() == B.cols())
{
if(block_size > std::min( N, std::min(M,K)) )
block_size = std::min( N, std::min(M,K));
std::vector<hsize_t> vsizetoRead, vstart,
vsizetoReadM, vstartM,
vsizetoReadK, vstartK;
getBlockPositionsSizes_hdf5( N, block_size, vstart, vsizetoRead );
getBlockPositionsSizes_hdf5( M, block_size, vstartM, vsizetoReadM );
getBlockPositionsSizes_hdf5( K, block_size, vstartK, vsizetoReadK );
// int ithreads = get_number_threads(threads, R_NilValue);
// int chunks = vstart.size()/ithreads;
#pragma omp parallel num_threads( get_number_threads(threads, R_NilValue) ) shared(A, B, C) //..// , chunk) private(tid )
{
#pragma omp for schedule (static) // collapse(3)
for (hsize_t ii = 0; ii < vstart.size(); ii ++)
{
for (hsize_t jj = 0; jj < vstartM.size(); jj++)
{
for (hsize_t kk = 0; kk < vstartK.size(); kk++)
{
C.block(vstart[ii], vstartM[jj], vsizetoRead[ii], vsizetoReadM[jj]) =
C.block(vstart[ii], vstartM[jj], vsizetoRead[ii], vsizetoReadM[jj]) +
( A.block(vstart[ii], vstartK[kk], vsizetoRead[ii], vsizetoReadK[kk]) *
B.block(vstartK[kk], vstartM[jj], vsizetoReadK[kk], vsizetoReadM[jj]) );
}
}
}
}
} else {
throw std::range_error("multiplication error: non-conformable arguments");
}
} catch(std::exception& ex) {
Rcpp::Rcerr<< "c++ exception Bblock_matrix_mul_parallel: "<<ex.what()<< " \n";
}
return(C);
}8 Usage Example
#include "BigDataStatMeth.hpp"
// Example usage
auto result = Bblock_matrix_mul_parallel(...);