hdf5_matrixVector_calculus
C++ Function Reference
1 Signature
BigDataStatMeth::hdf5Dataset * BigDataStatMeth::hdf5_matrixVector_calculus(BigDataStatMeth::hdf5Dataset *dsA, BigDataStatMeth::hdf5Dataset *dsB, BigDataStatMeth::hdf5Dataset *dsC, int function, bool bbyrows, bool bparal, Rcpp::Nullable< int > threads=R_NilValue)2 Description
Vector-matrix operations for HDF5 matrices.
3 Parameters
dsA(BigDataStatMeth::hdf5Dataset *): Input matrix datasetdsB(BigDataStatMeth::hdf5Dataset *): Input vector datasetdsC(BigDataStatMeth::hdf5Dataset *): Output matrix datasetfunction(int): Operation type (multiplication, addition, subtraction, division, power)bbyrows(bool): Whether to operate by rows or columnsbparal(bool): Whether to use parallel processingthreads(Rcpp::Nullable< int >): Number of threads for parallel processing
4 Returns
Type: BigDataStatMeth::hdf5Dataset *
5 Details
Performs vector-matrix operations on HDF5 datasets with support for parallel processing and row/column-wise operations.
6 Call Graph
7 Source Code
NoteImplementation
File: inst/include/hdf5Algebra/vectormatrix.hpp • Lines 208-335
inline BigDataStatMeth::hdf5Dataset* hdf5_matrixVector_calculus(
BigDataStatMeth::hdf5Dataset* dsA, BigDataStatMeth::hdf5Dataset* dsB,
BigDataStatMeth::hdf5Dataset* dsC, int function, bool bbyrows,
bool bparal, Rcpp::Nullable<int> threads = R_NilValue)
{
try{
std::vector<hsize_t> stride = {1, 1};
std::vector<hsize_t> block = {1, 1};
int blocksize = 0;
// Define blocksize atending number of elements in rows and cols
if( bbyrows == false) {
if( dsA->ncols() > MAXELEMSINBLOCK ) {
blocksize = 1;
} else {
hsize_t maxsize = std::max<hsize_t>( dsA->nrows(), dsA->ncols());
blocksize = std::ceil( MAXELEMSINBLOCK / maxsize);
}
} else {
if( dsA->nrows() > MAXELEMSINBLOCK) {
blocksize = 1;
} else {
hsize_t maxsize = std::max<hsize_t>( dsA->nrows(), dsA->ncols());
blocksize = std::ceil( MAXELEMSINBLOCK / maxsize);
}
}
dsC->createDataset( dsA->nrows(), dsA->ncols(), "real");
std::vector<double> vdB( dsB->nrows() * dsB->ncols());
dsB->readDatasetBlock( {0, 0}, { dsB->nrows(), dsB->ncols()}, stride, block, vdB.data() );
if( bbyrows == false) {
Eigen::RowVectorXd vWeights = Eigen::Map<Eigen::VectorXd, Eigen::Unaligned>(vdB.data(), vdB.size());
for(hsize_t i=0; (i * blocksize) <= dsA->nrows(); i++)
{
hsize_t sizetoread = 0;
if((i+1)*blocksize<dsA->nrows()) {
sizetoread = blocksize;
} else {
sizetoread = dsA->nrows()-(i*blocksize);
}
std::vector<hsize_t> offset = { i*blocksize, 0};
std::vector<hsize_t> count = {sizetoread, dsA->ncols()};
// Compute and write data
std::vector<double> vdA( count[0] * count[1]);
dsA->readDatasetBlock( {offset[0], offset[1]}, {count[0], count[1]}, stride, block, vdA.data() );
Eigen::Map<Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>> X (vdA.data(), count[0], count[1] );
if( function == 0 ) {
X = Rcpp_matrixVectorSum_byCol(X, vWeights);
} else if( function == 1 ) {
X = Rcpp_matrixVectorSubstract_byCol(X, vWeights);
} else if( function == 2 ) {
X = Rcpp_matrixVectorMultiplication_byCol(X, vWeights);
} else if( function == 3 ) {
X = Rcpp_matrixVectorDivision_byCol(X, vWeights);
} else if( function == 4 ) {
Rcpp_matrixVectorPower_byCol(X, vWeights);
}
dsC->writeDatasetBlock(Rcpp::wrap(X.transpose()), offset, count, stride, block, false);
}
} else {
Eigen::VectorXd vWeights = Eigen::Map<Eigen::VectorXd, Eigen::Unaligned>(vdB.data(), vdB.size());
for(hsize_t i=0; i*blocksize <= dsA->ncols() ; i++) {
hsize_t sizetoread = 0;
if( (i+1)*blocksize < dsA->ncols() ){
sizetoread = blocksize;
} else {
sizetoread = dsA->ncols()-(i*blocksize);
}
std::vector<hsize_t> offset = {0, i*blocksize};
std::vector<hsize_t> count = { dsA->nrows(), sizetoread};
// Compute and write data
std::vector<double> vdA( count[0] * count[1]);
dsA->readDatasetBlock( {offset[0], offset[1]}, {count[0], count[1]}, stride, block, vdA.data() );
Eigen::Map<Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor>> X (vdA.data(), count[0], count[1] );
if( function == 0 ) {
X = Rcpp_matrixVectorSum_byRow(X, vWeights);
} else if( function == 1 ) {
X = Rcpp_matrixVectorSubstract_byRow(X, vWeights);
} else if( function == 2 ) {
X = Rcpp_matrixVectorMultiplication_byRow(X, vWeights);
} else if( function == 3 ) {
X = Rcpp_matrixVectorDivision_byRow(X, vWeights);
} else if( function == 4 ) {
X = Rcpp_matrixVectorPow_byRow(X, vWeights);
}
offset = {i*blocksize, 0};
dsC->writeDatasetBlock(Rcpp::wrap(X.transpose()), offset, count, stride, block, true);
}
}
} catch( H5::FileIException& error ) { // catch failure caused by the H5File operations
// error.printErrorStack();
checkClose_file(dsA, dsB, dsC);
Rcpp::Rcerr<< "c++ exception hdf5_matrixVector_calculus (File IException)" ;
} catch( H5::DataSetIException& error ) { // catch failure caused by the DataSet operations
// error.printErrorStack();
checkClose_file(dsA, dsB, dsC);
Rcpp::Rcerr<<"c++ exception hdf5_matrixVector_calculus (DataSet IException)";
} catch(std::exception& error) {
checkClose_file(dsA, dsB, dsC);
Rcpp::Rcout<< "c++ exception vector-matrix functions: "<<error.what()<< " \n";
// return(dsC);
}
return(dsC);
}8 Usage Example
#include "BigDataStatMeth.hpp"
// Example usage
auto result = hdf5_matrixVector_calculus(...);