Doxygen 1.9.1
Toolkit for Adaptive Stochastic Modeling and Non-Intrusive ApproximatioN: Tasmanian v8.2 (development)
tsgLinearSolvers.hpp
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2017, Miroslav Stoyanov
3  *
4  * This file is part of
5  * Toolkit for Adaptive Stochastic Modeling And Non-Intrusive ApproximatioN: TASMANIAN
6  *
7  * Redistribution and use in source and binary forms, with or without modification, are permitted provided that the following conditions are met:
8  *
9  * 1. Redistributions of source code must retain the above copyright notice, this list of conditions and the following disclaimer.
10  *
11  * 2. Redistributions in binary form must reproduce the above copyright notice, this list of conditions
12  * and the following disclaimer in the documentation and/or other materials provided with the distribution.
13  *
14  * 3. Neither the name of the copyright holder nor the names of its contributors may be used to endorse
15  * or promote products derived from this software without specific prior written permission.
16  *
17  * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES,
18  * INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED.
19  * IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY,
20  * OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA,
21  * OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY,
22  * OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
23  *
24  * UT-BATTELLE, LLC AND THE UNITED STATES GOVERNMENT MAKE NO REPRESENTATIONS AND DISCLAIM ALL WARRANTIES, BOTH EXPRESSED AND IMPLIED.
25  * THERE ARE NO EXPRESS OR IMPLIED WARRANTIES OF MERCHANTABILITY OR FITNESS FOR A PARTICULAR PURPOSE, OR THAT THE USE OF THE SOFTWARE WILL NOT INFRINGE ANY PATENT,
26  * COPYRIGHT, TRADEMARK, OR OTHER PROPRIETARY RIGHTS, OR THAT THE SOFTWARE WILL ACCOMPLISH THE INTENDED RESULTS OR THAT THE SOFTWARE OR ITS USE WILL NOT RESULT IN INJURY OR DAMAGE.
27  * THE USER ASSUMES RESPONSIBILITY FOR ALL LIABILITIES, PENALTIES, FINES, CLAIMS, CAUSES OF ACTION, AND COSTS AND EXPENSES, CAUSED BY, RESULTING FROM OR ARISING OUT OF,
28  * IN WHOLE OR IN PART THE USE, STORAGE OR DISPOSAL OF THE SOFTWARE.
29  */
30 
31 #ifndef __TASMANIAN_LINEAR_SOLVERS_HPP
32 #define __TASMANIAN_LINEAR_SOLVERS_HPP
33 
35 
46 
62 namespace TasGrid{
63 
72 namespace TasmanianDenseSolver{
84  void solveLeastSquares(int n, int m, const double A[], double b[], double *x);
86  void solveLeastSquares(AccelerationContext const *acceleration, int n, int m, double A[], double b[], double *x);
93  template<typename scalar_type>
94  void solvesLeastSquares(AccelerationContext const *acceleration, int n, int m, scalar_type A[], int nrhs, scalar_type B[]);
96  template<typename scalar_type>
97  void solvesLeastSquaresGPU(AccelerationContext const *acceleration, int n, int m, scalar_type A[], int nrhs, scalar_type B[]);
98 }
99 
104 namespace TasmanianTridiagonalSolver{
106  std::vector<double> getSymmetricEigenvalues(int n, std::vector<double> const &diag, std::vector<double> const &offdiag);
107 
109  static constexpr int decompose_version = 1;
110 
115  void decompose(std::vector<double> &diag, std::vector<double> &off_diag, const double mu0, std::vector<double> &nodes,
116  std::vector<double> &weights);
117 
119  void decompose1(int n, std::vector<double> &d, std::vector<double> &e, std::vector<double> &z);
120 
122  void decompose2(std::vector<double> &diag, std::vector<double> &off_diag, const double mu0, std::vector<double> &nodes,
123  std::vector<double> &weights);
124 }
125 
129 namespace TasmanianFourierTransform{
133 
139  void fast_fourier_transform(std::vector<std::vector<std::complex<double>>> &data, std::vector<int> &num_points);
140 
144 
147  void fast_fourier_transform1D(std::vector<std::vector<std::complex<double>>> &data, std::vector<int> &indexes);
148 }
149 
153 namespace TasSparse{
154 
159 public:
161  WaveletBasisMatrix() : tol(Maths::num_tol), num_rows(0){}
163  WaveletBasisMatrix(AccelerationContext const *acceleration, const std::vector<int> &lpntr, const std::vector<std::vector<int>> &lindx, const std::vector<std::vector<double>> &lvals);
165  WaveletBasisMatrix(AccelerationContext const *acceleration, int cnum_rows, GpuVector<double> &&matrix)
166  : tol(Maths::num_tol), num_rows(cnum_rows), gpu_dense(std::move(matrix)) { factorize(acceleration); }
168  ~WaveletBasisMatrix() = default;
169 
178 
180  static bool useDense(AccelerationContext const *acceleration, int nrows){
182  and not (acceleration->algorithm_select == AccelerationContext::algorithm_autoselect and nrows > 10000));
183  }
184 
186  bool isSparse() const{ return (num_rows > 0) and dense.empty(); }
188  bool isDense() const{ return (num_rows > 0) and not dense.empty(); }
189 
191  int getNumRows() const{ return num_rows; }
192 
194  void invertTransposed(AccelerationContext const *acceleration, double b[]) const;
195 
197  void invert(AccelerationContext const *acceleration, int num_colums, double B[]);
198 
200  template<bool transpose, bool blas> void solve(const double b[], double x[]) const;
201 
202 protected:
204  void factorize(AccelerationContext const *acceleration);
206  void computeILU();
208  template<bool transpose> void applyILU(double x[]) const;
210  template<bool transpose> void apply(double const x[], double r[]) const;
212  void residual(double const x[], double const b[], double r[]) const;
214  static constexpr bool use_transpose = true;
216  static constexpr bool no_transpose = false;
218  static constexpr bool use_blas = true;
220  static constexpr bool no_blas = false;
221 
222 private:
223  double tol;
224  int num_rows;
225  std::vector<int> pntr, indx, indxD;
226  std::vector<double> vals, ilu;
227  std::vector<double> dense; // if using the dense format
228  std::vector<int> ipiv; // pivots for the dense factorize
229 
230  GpuVector<double> gpu_dense; // GPU factors
231  GpuVector<int_gpu_lapack> gpu_ipiv;
232 };
233 
234 }
235 
236 }
237 
238 #endif
Used to manipulate the wavelet values and solve for the wavelet coefficients.
Definition: tsgLinearSolvers.hpp:158
WaveletBasisMatrix()
Default constructor, create an empty matrix.
Definition: tsgLinearSolvers.hpp:161
void invert(AccelerationContext const *acceleration, int num_colums, double B[])
Overwrites the row-major matrix B with X that solves .
WaveletBasisMatrix(AccelerationContext const *acceleration, const std::vector< int > &lpntr, const std::vector< std::vector< int >> &lindx, const std::vector< std::vector< double >> &lvals)
Initialize the matrix with the given set of indexes.
bool isDense() const
Return true if using the dense mode and false is using sparse mode or empty.
Definition: tsgLinearSolvers.hpp:188
void factorize(AccelerationContext const *acceleration)
Computes the factorization of the matrix, ILU or PLU using the provided acceleration context (called ...
WaveletBasisMatrix(WaveletBasisMatrix &&)=default
Move constructor.
static constexpr bool no_blas
Expressive call to no-blas variant method.
Definition: tsgLinearSolvers.hpp:220
void solve(const double b[], double x[]) const
Solve op(A) x = b where op is either identity (find the coefficients) or transpose (find the interpol...
WaveletBasisMatrix(WaveletBasisMatrix const &)=delete
Cannot copy.
void computeILU()
Compute the incomplete lower-upper decomposition of the matrix (zero extra fill).
WaveletBasisMatrix(AccelerationContext const *acceleration, int cnum_rows, GpuVector< double > &&matrix)
Initialize the matrix in dense mode with the given data.
Definition: tsgLinearSolvers.hpp:165
void invertTransposed(AccelerationContext const *acceleration, double b[]) const
Overwrites vector b with x that solves .
bool isSparse() const
Return true if using the sparse mode and false is using dense mode or empty.
Definition: tsgLinearSolvers.hpp:186
void residual(double const x[], double const b[], double r[]) const
Computes the residual.
int getNumRows() const
Return the number of rows in the matrix.
Definition: tsgLinearSolvers.hpp:191
static constexpr bool no_transpose
Expressive call to no transpose method.
Definition: tsgLinearSolvers.hpp:216
void applyILU(double x[]) const
Apply the preconditioner on a vector.
static constexpr bool use_blas
Expressive call to blas variant method.
Definition: tsgLinearSolvers.hpp:218
void apply(double const x[], double r[]) const
Sets to be the product of the sparse matrix times x.
static constexpr bool use_transpose
Expressive call to transpose method.
Definition: tsgLinearSolvers.hpp:214
static bool useDense(AccelerationContext const *acceleration, int nrows)
Decide between sparse and dense variant of the algorithms.
Definition: tsgLinearSolvers.hpp:180
WaveletBasisMatrix & operator=(WaveletBasisMatrix const &)=delete
Cannot copy.
~WaveletBasisMatrix()=default
Default destructor.
void fast_fourier_transform(std::vector< std::vector< std::complex< double >>> &data, std::vector< int > &num_points)
Transfrom the data for a multi-dimensional tensor.
void solveLeastSquares(int n, int m, const double A[], double b[], double *x)
Least squares solver, used to infer anisotropic coefficients and thus rarely exceeds 10 - 20.
void fast_fourier_transform1D(std::vector< std::vector< std::complex< double >>> &data, std::vector< int > &indexes)
Perform one dimensional fast-fourier-transform (using radix-3).
constexpr double num_tol
Numerical tolerance for various algorithms.
Definition: tsgMathUtils.hpp:132
void solvesLeastSquares(AccelerationContext const *acceleration, int n, int m, scalar_type A[], int nrhs, scalar_type B[])
Least squares solver, operates on multiple right-hand sides and row-major matrices.
void solvesLeastSquaresGPU(AccelerationContext const *acceleration, int n, int m, scalar_type A[], int nrhs, scalar_type B[])
Overload that accepts arrays on the GPU device.
void decompose(std::vector< double > &diag, std::vector< double > &off_diag, const double mu0, std::vector< double > &nodes, std::vector< double > &weights)
Method for tridiagonal eigenvalue decomposition, used to compute nodes and weights for Gaussian rules...
void decompose1(int n, std::vector< double > &d, std::vector< double > &e, std::vector< double > &z)
(decompose_version == 1) TASMANIAN's first internal implementation.
std::vector< double > getSymmetricEigenvalues(int n, std::vector< double > const &diag, std::vector< double > const &offdiag)
Method for computing the eigenvalues of a symmetric matrix in place, using an LAPACK wrapper.
void decompose2(std::vector< double > &diag, std::vector< double > &off_diag, const double mu0, std::vector< double > &nodes, std::vector< double > &weights)
(decompose_version == 2) Based on the ALGOL code for Golub's 1967 report "Calculation of Gauss Quadra...
Encapsulates the Tasmanian Sparse Grid module.
Definition: TasmanianSparseGrid.hpp:68
Wrapper class around GPU device ID, acceleration type and GpuEngine.
Definition: tsgAcceleratedDataStructures.hpp:576
AlgorithmPreference algorithm_select
The preference to use dense or sparse algorithms.
Definition: tsgAcceleratedDataStructures.hpp:609
@ algorithm_sparse
Use sparse algorithm.
Definition: tsgAcceleratedDataStructures.hpp:582
@ algorithm_autoselect
Use automatically select based on heuristics.
Definition: tsgAcceleratedDataStructures.hpp:584
Data structures for interacting with CUDA and MAGMA environments.