This is the main HODLR class that includes all the fast HODLR solvers. More...
#include <HODLR_Matrix.hpp>
Public Member Functions | |
HODLR_Matrix (Eigen::MatrixXd &inputMatrix, int inputSizeThreshold=30, std::string LR_Method="partialPiv_ACA") | |
This constructor initializes the class with a dense matrix and an optional leaf size threshold. More... | |
HODLR_Matrix (Eigen::MatrixXd &inputMatrix, int inputSizeThreshold, user_IndexTree &input_IndexTree) | |
This constructor initializes the class with a dense matrix and a user specified indexing scheme which will be used to create the HODLR index tree. This constructor initializes the class with a dense matrix. More... | |
HODLR_Matrix (Eigen::MatrixXd &inputMatrix, Eigen::SparseMatrix< double > &inputGraph, int inputSizeThreshold=30) | |
This constructor initializes the class with a dense matrix and an interaction graph (sparse matrix). The grpah is going to be used primarily as the interaction graph in the BDLR scheme. More... | |
HODLR_Matrix (Eigen::MatrixXd &inputMatrix, Eigen::SparseMatrix< double > &inputGraph, int inputSizeThreshold, user_IndexTree &input_IndexTree) | |
This constructor initializes the class with a dense matrix and an interaction graph (sparse matrix) and a user defined indexing schemes. The grpah is going to be used primarily as the interaction graph in the BDLR scheme. The user defined indexing scheme will be used to create the HODLR index tree. More... | |
HODLR_Matrix (Eigen::SparseMatrix< double > &inputMatrix, int inputSizeThreshold=30, std::string LR_Method="PS_Sparse") | |
This constructor initializes the class with a sparse matrix and optional leaf size and low-rank approximation method parameters. More... | |
HODLR_Matrix (Eigen::SparseMatrix< double > &inputMatrix, int inputSizeThreshold, user_IndexTree &input_IndexTree, std::string LR_Method="PS_Sparse") | |
This constructor initializes the class with a sparse matrix and a user defined indexing scheme used to create the the HODLR indexing tree. More... | |
HODLR_Matrix (Eigen::SparseMatrix< double > &inputMatrix, Eigen::SparseMatrix< double > &inputGraph, int inputSizeThreshold, user_IndexTree &input_IndexTree, std::string LR_Method="PS_Sparse") | |
This constructor initializes the class with a sparse matrix and an interaction graph (sparse matrix) and a user defined indexing schemes. The grpah is going to be used primarily as the interaction graph in the BDLR scheme. The user defined indexing scheme will be used to create the HODLR index tree. More... | |
HODLR_Matrix (int numRows, int numCols, double(*inputKernel)(int i, int j, void *inputKernelData), void *inputKernelData, int inputSizeThreshold) | |
HODLR_Matrix (int numRows, int numCols, double(*inputKernel)(int i, int j, void *inputKernelData), void *inputKernelData, Eigen::SparseMatrix< double > &inputGraph, int inputSizeThreshold) | |
HODLR_Matrix (int numRows, int numCols, double(*inputKernel)(int i, int j, void *inputKernelData), void *inputKernelData, int inputSizeThreshold, user_IndexTree &input_IndexTree) | |
HODLR_Matrix (int numRows, int numCols, double(*inputKernel)(int i, int j, void *inputKernelData), void *inputKernelData, Eigen::SparseMatrix< double > &inputGraph, int inputSizeThreshold, user_IndexTree &input_IndexTree) | |
HODLR_Matrix (const HODLR_Matrix &rhs) | |
void | storeLRinTree () |
Eigen::MatrixXd | recLU_Solve (const Eigen::MatrixXd &input_RHS) |
Eigen::MatrixXd | recSM_Solve (const Eigen::MatrixXd &input_RHS) |
void | recLU_Compute () |
Eigen::MatrixXd | extendedSp_Solve (const Eigen::MatrixXd &input_RHS) |
Eigen::MatrixXd | iterative_Solve (const Eigen::MatrixXd &input_RHS, const int maxIterations, const double stop_tolerance, const double init_LRTolerance, const std::string input_LR_Method, const std::string directSolve_Method) |
void | set_LRTolerance (double tolerance) |
void | set_minPivot (double minPivot) |
void | set_LRMethod (std::string input_LRMethod) |
void | set_FreeMatrixMemory (bool inputVal) |
void | set_BoundaryDepth (int inputBoundaryDepth) |
void | set_recLUFactorizedFlag (bool factorized) |
void | set_LeafConst () |
void | saveExtendedSp (std::string savePath) |
double | get_recLU_FactorizationTime () const |
double | get_recLU_SolveTime () const |
double | get_recLU_TotalTime () const |
double | get_extendedSp_AssemblyTime () const |
double | get_extendedSp_FactorizationTime () const |
double | get_extendedSp_SolveTime () const |
double | get_extendedSp_TotalTime () const |
double | get_LR_ComputationTime () const |
double | get_totalIter_SolveTime () const |
double | get_MatrixSize () const |
int | rows () const |
int | cols () const |
double | norm () |
double | determinant () |
double | logAbsDeterminant () |
Eigen::MatrixXd | block (int min_i, int min_j, int numRows, int numCols) |
Eigen::MatrixXd | row (int row) |
Eigen::MatrixXd | col (int col) |
HODLR_Matrix | topDiag () |
HODLR_Matrix | bottDiag () |
void | keepTopDiag () |
void | keepBottDiag () |
Eigen::MatrixXd & | returnTopOffDiagU () |
Eigen::MatrixXd & | returnTopOffDiagV () |
Eigen::MatrixXd & | returnTopOffDiagK () |
Eigen::MatrixXd & | returnBottOffDiagU () |
Eigen::MatrixXd & | returnBottOffDiagV () |
Eigen::MatrixXd & | returnBottOffDiagK () |
void | check_Structure () |
double | calcAbsDiff () |
Eigen::MatrixXd | createExactHODLR (const int rank, int input_MatrixSize, const int inpt_SizeThreshold) |
void | saveSolverInfo (const std::string outputFileName) |
void | freeMatrixData () |
void | destroyAllData () |
void | recalculateSize () |
void | correctIndices () |
void | initInfoVectors () |
Public Attributes | |
bool | printLevelRankInfo |
bool | printLevelAccuracy |
bool | printLevelInfo |
bool | printResultInfo |
Friends | |
void | splitAtTop (HODLR_Matrix &self, HODLR_Matrix &topHODLR, HODLR_Matrix &bottHODLR) |
void | extendAddUpdate (HODLR_Matrix &parentHODLR, std::vector< Eigen::MatrixXd * > D_Array, std::vector< HODLR_Matrix * > D_HODLR_Array, std::vector< Eigen::MatrixXd * > U_Array, std::vector< Eigen::MatrixXd * > V_Array, std::vector< std::vector< int > > &updateIdxVec_Array_D, std::vector< std::vector< int > > &updateIdxVec_Array_D_HODLR, double tol, std::string mode) |
void | extendAddUpdate (HODLR_Matrix &parentHODLR, HODLR_Matrix &D_HODLR, std::vector< int > &updateIdxVec, double tol, std::string mode) |
void | extendAddUpdate (HODLR_Matrix &parentHODLR, Eigen::MatrixXd &U, Eigen::MatrixXd &V, std::vector< int > &updateIdxVec, double tol, std::string mode) |
HODLR_Matrix | extend (std::vector< int > &extendIdxVec, int parentSize, HODLR_Matrix &childHODLR) |
void | extendAddUpdate (HODLR_Matrix &parentHODLR, Eigen::MatrixXd &D, std::vector< int > &updateIdxVec, double tol, std::string mode) |
This is the main HODLR class that includes all the fast HODLR solvers.
This class can be used as a sibling of the Eigen::MatrixXd class in cases where access to matrix entries is required.