Rcpp_InvCholesky_hdf5
C++ Function Reference
1 Signature
void BigDataStatMeth::Rcpp_InvCholesky_hdf5(BigDataStatMeth::hdf5Dataset *inDataset, BigDataStatMeth::hdf5DatasetInternal *outDataset, bool bfull, long dElementsBlock, Rcpp::Nullable< int > threads)2 Description
Computes matrix inverse using Cholesky decomposition with HDF5 storage.
3 Parameters
inDataset(BigDataStatMeth::hdf5Dataset *): Input matrix dataset (must be symmetric positive-definite)outDataset(BigDataStatMeth::hdf5DatasetInternal *): Output dataset for the inverse matrixbfull(bool): If true, computes full matrix inverse; if false, only lower triangular partdElementsBlock(long): Block size for processing (minimum 2 * matrix dimension)threads(Rcpp::Nullable< int >): Number of threads for parallel processing (optional)
4 Details
inDatasetInput matrix dataset (must be symmetric positive-definite) outDatasetOutput dataset for the inverse matrix bfullIf true, computes full matrix inverse; if false, only lower triangular part dElementsBlockBlock size for processing (minimum 2 * matrix dimension) threadsNumber of threads for parallel processing (optional) This function performs matrix inversion in three steps:Cholesky decomposition (A = LL^T)Inverse of the Cholesky factor (L^-1)Computation of full inverse (A^-1 = L^-T L^-1)
5 Call Graph
6 Source Code
NoteImplementation
File: inst/include/hdf5Algebra/matrixInvCholesky.hpp • Lines 260-315
inline void Rcpp_InvCholesky_hdf5 ( BigDataStatMeth::hdf5Dataset* inDataset,
BigDataStatMeth::hdf5DatasetInternal* outDataset,
bool bfull, long dElementsBlock,
Rcpp::Nullable<int> threads = R_NilValue )
{
try{
int nrows = inDataset->nrows();
int ncols = inDataset->ncols();
// Rcpp::Rcout<<"\nIniciem cholesky";
int res = Cholesky_decomposition_hdf5(inDataset, outDataset, nrows, ncols, dElementsBlock, threads);
if(res == 0)
{
// Rcpp::Rcout<<"\nDins res==0 -> anem Inverse_of_Cholesky_decomposition_hdf5...";
Inverse_of_Cholesky_decomposition_hdf5( outDataset, nrows, ncols, dElementsBlock, threads);
// Rcpp::Rcout<<"\nDins res==0 -> anem Inverse_Matrix_Cholesky_parallel...";
Inverse_Matrix_Cholesky_hdf5( outDataset, nrows, ncols, dElementsBlock, threads);
if( bfull==true ) {
setUpperTriangularMatrix( outDataset, dElementsBlock);
}
}
} catch( H5::FileIException& error ) {
checkClose_file(inDataset, outDataset);
inDataset = outDataset = nullptr;
Rcpp::Rcerr<<"c++ exception Rcpp_InvCholesky_hdf5 (File IException)";
return void();
} catch( H5::GroupIException & error ) {
checkClose_file(inDataset, outDataset);
inDataset = outDataset = nullptr;
Rcpp::Rcerr << "c++ exception Rcpp_InvCholesky_hdf5 (Group IException)";
return void();
} catch( H5::DataSetIException& error ) {
checkClose_file(inDataset, outDataset);
inDataset = outDataset = nullptr;
Rcpp::Rcerr << "c++ exception Rcpp_InvCholesky_hdf5 (DataSet IException)";
return void();
} catch(std::exception& ex) {
checkClose_file(inDataset, outDataset);
inDataset = outDataset = nullptr;
Rcpp::Rcerr << "c++ exception Rcpp_InvCholesky_hdf5" << ex.what();
return void();
} catch (...) {
checkClose_file(inDataset, outDataset);
inDataset = outDataset = nullptr;
Rcpp::Rcerr<<"\nC++ exception Rcpp_InvCholesky_hdf5 (unknown reason)";
return void();
}
return void();
}7 Usage Example
#include "BigDataStatMeth.hpp"
// Example usage
auto result = Rcpp_InvCholesky_hdf5(...);