Rcpp_vector_multiply_hdf5
C++ Function Reference
1 Signature
BigDataStatMeth::hdf5Dataset * BigDataStatMeth::Rcpp_vector_multiply_hdf5(BigDataStatMeth::hdf5Dataset *dsA, BigDataStatMeth::hdf5Dataset *dsB, BigDataStatMeth::hdf5Dataset *dsC, bool bparal, Rcpp::Nullable< int > threads=R_NilValue)2 Description
Pure vector multiplication for HDF5 vectors.
3 Parameters
dsA(BigDataStatMeth::hdf5Dataset *): First input vector datasetdsB(BigDataStatMeth::hdf5Dataset *): Second input vector datasetdsC(BigDataStatMeth::hdf5Dataset *): Output vector datasetbparal(bool): Whether to use parallel processingthreads(Rcpp::Nullable< int >): Number of threads for parallel processing (optional)
4 Returns
Pointer to result dataset
5 Details
Performs optimized element-wise multiplication C = A * B where A, B, and C are HDF5 vector datasets. Uses same optimization strategy as vector addition.
6 Call Graph
7 Source Code
NoteImplementation
File: inst/include/hdf5Algebra/vectorOperations.hpp • Lines 263-323
inline BigDataStatMeth::hdf5Dataset* Rcpp_vector_multiply_hdf5(
BigDataStatMeth::hdf5Dataset* dsA, BigDataStatMeth::hdf5Dataset* dsB, BigDataStatMeth::hdf5Dataset* dsC,
bool bparal, Rcpp::Nullable<int> threads = R_NilValue)
{
try {
hsize_t sizeA = validateVector(dsA);
hsize_t sizeB = validateVector(dsB);
if (sizeA == 0 || sizeB == 0) {
Rcpp::Rcout << "vector multiply error: inputs are not vectors\n";
return dsC;
}
if (sizeA != sizeB) {
Rcpp::Rcout << "vector multiply error: non-conformable vector dimensions\n";
return dsC;
}
hsize_t rowsA = dsA->nrows();
hsize_t colsA = dsA->ncols();
dsC->createDataset(colsA, rowsA, "real");
std::vector<hsize_t> stride = {1, 1};
std::vector<hsize_t> block = {1, 1};
std::vector<double> vdA(sizeA);
std::vector<double> vdB(sizeB);
dsA->readDatasetBlock({0, 0}, {rowsA, colsA}, stride, block, vdA.data());
dsB->readDatasetBlock({0, 0}, {dsB->nrows(), dsB->ncols()}, stride, block, vdB.data());
if (bparal && sizeA > 10000) {
#pragma omp parallel num_threads(get_threads(bparal, threads))
{
#pragma omp for schedule(static)
for (hsize_t i = 0; i < sizeA; ++i) {
vdA[i] *= vdB[i];
}
}
} else {
std::transform(vdA.begin(), vdA.end(), vdB.begin(), vdA.begin(), std::multiplies<double>());
}
dsC->writeDatasetBlock(vdA, {0, 0}, {rowsA, colsA}, stride, block);
} catch(H5::FileIException& error) {
checkClose_file(dsA, dsB, dsC);
Rcpp::Rcerr << "\nc++ exception Rcpp_vector_multiply_hdf5 (File IException)";
// return dsC;
} catch(H5::DataSetIException& error) {
checkClose_file(dsA, dsB, dsC);
Rcpp::Rcerr << "\nc++ exception Rcpp_vector_multiply_hdf5 (DataSet IException)";
// return dsC;
} catch(std::exception& ex) {
checkClose_file(dsA, dsB, dsC);
Rcpp::Rcerr << "\nc++ exception Rcpp_vector_multiply_hdf5: " << ex.what();
// return dsC;
}
return dsC;
}8 Usage Example
#include "BigDataStatMeth.hpp"
// Example usage
auto result = Rcpp_vector_multiply_hdf5(...);