HOLDER_PAC
 All Classes Functions Modules
helperFunctions.hpp
1 #ifndef HELPERFUNCTIONS_HODLR_SOLVER_HPP
2 #define HELPERFUNCTIONS_HODLR_SOLVER_HPP
3 
4 //C++ Dependencies
5 #include <algorithm>
6 #include <chrono>
7 #include <cmath>
8 #include <ctime>
9 #include <fstream>
10 #include <iostream>
11 #include <random>
12 #include <vector>
13 
14 //External Dependencies
15 #include <Eigen/Dense>
16 #include <Eigen/Sparse>
17 
18 //Custom Dependencies
19 #include "HODLR_Matrix.hpp"
20 #include "user_IndexTree.hpp"
21 
22 /******************Pre-programmed radial basis function kernels:**********************/
23 
24 //1+r^2
25 double quadraticKernel(const double r);
26 
27 //sqrt(1+r^2)
28 double multiQuadraticKernel(const double r);
29 
30 //1/(1+r^2)
31 double inverseQuadraticKernel(const double r);
32 
33 //1/sqrt(1+r^2)
34 double inverseMultiQuadraticKernel(const double r);
35 
36 //exp(-r^2)
37 double gaussianKernel(const double r);
38 
39 //exp(-r)
40 double exponentialKernel(const double r);
41 
42 //log(1+r)
43 double logarithmicKernel(const double r);
44 
45 //1/r
46 double oneOverRKernel(const double r);
47 
48 //1/sqrt(r)
49 double oneOverSqRKernel(const double r);
50 
51 //log(r)
52 double logRKernel(const double r);
53 
54 /**************************************************************************************/
55 
56 
57 std::vector<int> createUniqueRndIdx(const int min, const int max,const int n);
58 
59 std::vector<int> createSequentialVec(const int min,const int size);
60 
61 /* Function: makeMatrixFrom1DInterval
62  * ----------------------------------
63  * This function creates a dense interaction matrix given a set of row and column points given a kernel.
64  * The points must lie on a 1D manifold.
65  * rowPts : Coordinates of the row points.
66  * colPts : Coordinates of the column points.
67  * diagValue : The diagonal entry value of the dense matrix.
68  * kernel : Pointer to the kernel function.
69  */
70 Eigen::MatrixXd makeMatrixFrom1DInterval(const std::vector<double> rowPts,const std::vector<double> colPts,const double diagValue ,double (*kernel)(const double));
71 
72 
73 /* Function: makeMatrix1DUniformPts
74  * --------------------------------
75  * This function creates a dense interaction matrix arising from the interaction of uniform points on a 1D interval.
76  * This function is mostly a wrapper for makeMatrixFrom1DInterval.
77  * minRowPt : Lower bound of the 1D interval for the row points.
78  * maxRowPt : Upper bound for the 1D interval for the row points.
79  * minColPt : Lower bound of the 1D interval for the column points.
80  * maxColpt : Upper bound of the 1D interval for the row points.
81  * numRows : Number of rows (row points).
82  * numCols : Number of columns (column points).
83  * diagValue : The diagonal entry value of the dense matrix.
84 `* kernel : Pointer to the kernel function.
85  */
86 Eigen::MatrixXd makeMatrix1DUniformPts(const int minRowPt, const int maxRowPt, const int minColPt, const int maxColPt, const int numRows, const int numCols, const double diagValue, double (*kernel)(const double));
87 
88 
89 /* Function: testACASolverConv1DUniformPts
90  * ---------------------------------------
91  * This function creates a convergence plot of solver relative error vs ACA tolerance for a dense self interaction matrix.
92  * The test dense matrix, is an interaction matrix arising from the self interaction of uniform points on a 1D interval.
93  * intervalMin : Lower bound of the 1D interval.
94  * intervalMax : Upper bound of the 1D interval.
95  * numPts : Number of interacting points (matrix size).
96  * diagValue : The diagonal entry value of the dense matrix.
97  * exactSoln : The test right hand side of the linear system.
98  * outputFileName : Path of the output file.
99  * kernel : Pointer to the kernel function.
100  * solverType : Type of HODLR solver. Default is recLU.
101  */
102 void testACASolverConv1DUnifromPts(const double intervalMin,const double intervalMax, const int numPts, const double diagValue, Eigen::VectorXd exactSoln, std::string outputFileName, double (*kernel)(const double r), std::string solverType = "recLU");
103 
104 
105 /* Function: testACASolverSpeed1DUniformPts
106  * ---------------------------------------
107  * This function creates a plot of solver CPU Time vs matrix size for dense self interaction matrice with various sizes.
108  * The test dense matrices, are interaction matrix arising from the self interaction of uniform points on a 1D interval.
109  * intervalMin : Lower bound of the 1D interval.
110  * intervalMax : Upper bound of the 1D interval.
111  * diagValue : The diagonal entry value of the dense matrix.
112  * ACATolerance : The ACA tolerance of the HODLR solver.
113  * outputFileName : Path of the output file.
114  * kernel : Pointer to the kernel function.
115  * solverType : Type of HODLR solver. Default is recLU.
116  */
117 void testACASolverSpeed1DUniformPts(const double intervalMin, const double intervalMax,const double diagValue,const double ACATolerance, std::string outputFileName, double (*kernel)(const double), std::string solverType = "recLU");
118 
119  /* Function: testACASolverSpeed1DUniformPts_FixedSize
120  * ---------------------------------------
121  * This function ccalculates the cpu time of the various stages in a solve process for a fixed-size matrix.
122  * The test dense matrices, are interaction matrix arising from the self interaction of uniform points on a 1D interval.
123  * intervalMin : Lower bound of the 1D interval.
124  * intervalMax : Upper bound of the 1D interval.
125  * diagValue : The diagonal entry value of the dense matrix.
126  * LR_Tolerance : The ACA tolerance of the HODLR solver.
127  * outputFileName : Path of the output file.
128  * kernel : Pointer to the kernel function.
129  * matrixSize : Size of the matrix.
130  * solverType : Type of HODLR solver. Default is recLU.
131  */
132 void testACASolverSpeed1DUniformPts_FixedSize(const double intervalMin, const double intervalMax,const double diagValue,const double LR_Tolerance, std::string outputFileName, double (*kernel)(const double), const int matrixSize, std::string solverType);
133 
134 double eigenPartialPivLUSpeed(const Eigen::MatrixXd & inputMatrix);
135 
136 void testSolverSpeed(const std::string inputFilePath,const std::string outputFilePath,const int sizeThreshold,std::string solverType,user_IndexTree & usrTree);
137 
138 void testBoundaryLRSolver(const std::string inputMatrixFileName,const std::string inputGraphFileName,const std::string outputFileName,const double iterInitTol,const int sizeThreshold,const int depth);
139 
140 void testBoundaryLRSolver(const std::string inputMatrixFileName,const std::string inputGraphFileName,const std::string outputFileName,const double iterInitTol,const int sizeThreshold,const int depth,user_IndexTree &usrTree);
141 
142 void analyzeRank(const std::string inputMatrixFileName,const std::string inputGraphFileName,const std::string outputFileName,const int input_Min_i,const int input_Min_j,const int input_NumRows, const int input_NumCols,std::string mode = "default");
143 
144 #endif
Definition: user_IndexTree.hpp:7