prepareForParallelization

C++ Function Reference

1 Signature

std::vector< svdPositions > BigDataStatMeth::prepareForParallelization(T *dsA, int M, int k, bool transp, int block_size, std::string strDatasetName)

2 Parameters

  • dsA (T *)
  • M (int)
  • k (int)
  • transp (bool)
  • block_size (int)
  • strDatasetName (std::string)

3 Returns

Type: class T

4 Call Graph

Function dependencies

5 Source Code

File: inst/include/hdf5Algebra/matrixSvdBlock.hppLines 136-230

std::vector<svdPositions> prepareForParallelization( T* dsA, int M, int k, bool transp, int block_size, std::string strDatasetName)
{

    static_assert(std::is_same<T*, BigDataStatMeth::hdf5Dataset* >::value ||
                  std::is_same<T*, BigDataStatMeth::hdf5DatasetInternal* >::value,
                  "Error - type not allowed");


    std::vector<svdPositions> pos;
    BigDataStatMeth::hdf5Dataset* unlimDataset = nullptr;

    // Rcpp::Rcout<<"\nAnem a preparar per paralelitzar... a veure que fem per aquí perquè no m'agrada massa... \n";

    try{

        // int realsizeread, cummoffset;
        int realsizeread;
        int irows = dsA->ncols();
        int icols = dsA->nrows();

        
        if(transp == false) {
            irows = dsA->ncols();
            icols = dsA->nrows();    
        } else {
            irows = dsA->nrows();
            icols = dsA->ncols();    
        }
        

        for( int i = 0; i< M ; i++)
        {
            int maxsizetoread = block_size;

            pos.push_back(svdPositions());

            pos[i].strDatasetName = strDatasetName + std::to_string(i/(M/k));
            pos[i].totOffset = getInitialPosition( transp, (unsigned long long)(i*block_size) ); // Initial read position

            // Get max block size to read - for blocks smaller than default block size
            if( ((i+1)*block_size) > icols)
                maxsizetoread = icols - (i*block_size);

            if( i+1 == M && icols - maxsizetoread!=0) {
                realsizeread = icols - (i*block_size);
            } else{
                realsizeread = maxsizetoread;
            }

            pos[i].count = getSizetoRead(transp, (unsigned long long)(realsizeread), icols, irows );

            if( i%(M/k) == 0 || ( (i%(M/k) > 0 &&  !exists_HDF5_element(dsA->getFileptr(),  pos[i].strDatasetName)) ) )
            {
                pos[i].cummoffset = 0;
            } else {
                pos[i].cummoffset = pos[i-1].cummoffset + pos[i-1].count[0];
            }

            pos[i].partOffset[1] = pos[i].cummoffset;
        }


    } catch( H5::FileIException& error ) {
        checkClose_file(dsA, unlimDataset);
        Rcpp::Rcerr<<"\nc++ exception prepareForParallelization (File IException)\n";
        return(pos);
    } catch( H5::DataSetIException& error ) {
        checkClose_file(dsA, unlimDataset);
        Rcpp::Rcerr<<"\nc++ exception prepareForParallelization (DataSet IException)\n";
        return(pos);
    } catch( H5::GroupIException& error ) {
        checkClose_file(dsA, unlimDataset);
        Rcpp::Rcerr<<"\nc++ exception prepareForParallelization (Group IException)\n";
        return(pos);
    } catch( H5::DataTypeIException& error ) {
        checkClose_file(dsA, unlimDataset);
        Rcpp::Rcerr<<"\nc++ exception prepareForParallelization (DataType IException)\n";
        return(pos);
    } catch( H5::DataSpaceIException& error ) {
        checkClose_file(dsA, unlimDataset);
        Rcpp::Rcerr<<"\nc++ exception prepareForParallelization (DataSpace IException)\n";
        return(pos);
    } catch(std::exception &ex) {
        checkClose_file(dsA, unlimDataset);
        Rcpp::Rcerr<<"c++ exception prepareForParallelization \n"<< ex.what();
        return(pos);
    } catch (...) {
        checkClose_file(dsA, unlimDataset);
        Rcpp::Rcerr<<"\nC++ exception prepareForParallelization (unknown reason)";
        return(pos);
    }


    return(pos);
}

6 Usage Example

#include "BigDataStatMeth.hpp"

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