TBTK
Need a break? Support the development by playing Polarity Puzzles
LUSolver.h
Go to the documentation of this file.
1 /* Copyright 2017 Kristofer Björnson
2  *
3  * Licensed under the Apache License, Version 2.0 (the "License");
4  * you may not use this file except in compliance with the License.
5  * You may obtain a copy of the License at
6  *
7  * http://www.apache.org/licenses/LICENSE-2.0
8  *
9  * Unless required by applicable law or agreed to in writing, software
10  * distributed under the License is distributed on an "AS IS" BASIS,
11  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
12  * See the License for the specific language governing permissions and
13  * limitations under the License.
14  */
15 
23 #ifndef COM_DAFER45_TBTK_LU_SOLVER
24 #define COM_DAFER45_TBTK_LU_SOLVER
25 
26 #include "TBTK/Communicator.h"
27 #include "TBTK/Matrix.h"
28 #include "TBTK/SparseMatrix.h"
29 
30 #include <complex>
31 
32 #include "slu_zdefs.h"
33 
34 namespace TBTK{
35 
37 class LUSolver : public Communicator{
38 public:
42  enum class DataType {None, Double, ComplexDouble};
43 
45  LUSolver();
46 
48  ~LUSolver();
49 
54  void setMatrix(const SparseMatrix<double> &sparseMatrix);
55 
60  void setMatrix(const SparseMatrix<std::complex<double>> &sparseMatrix);
61 
69 
74  void solve(Matrix<double> &b);
75 
80  void solve(Matrix<std::complex<double>> &b);
81 private:
83  SuperMatrix *L;
84 
86  SuperMatrix *U;
87 
89  int *rowPermutations;
90 
92  int *columnPermutations;
93 
95  SuperLUStat_t *statistics;
96 
98  DataType matrixDataType;
99 
101  void allocatePermutationMatrices(
102  unsigned int numRows,
103  unsigned int numColumns
104  );
105 
107  void initStatistics();
108 
110  void allocateLUMatrices();
111 
113  void initOptionsAndPermutationMatrices(
114  superlu_options_t &options,
115  SuperMatrix &matrix
116  );
117 
119  void performLUFactorization(SuperMatrix &matrix);
120 
122  void checkSolveAssert(unsigned int numRows);
123 
125  void checkXgstrfErrors(
126  int info,
127  std::string functionName,
128  int numColumns
129  );
130 
132  void checkXgstrsErrors(int info, std::string functionName);
133 };
134 
136  return matrixDataType;
137 }
138 
139 }; //End of namespace TBTK
140 
141 #endif
void setMatrix(const SparseMatrix< double > &sparseMatrix)
Sparse matrix.
Custom matrix.
Definition: Matrix.h:33
Definition: Boolean.h:32
Definition: SparseMatrix.h:35
Base class for classes that can communicate their status during execution.
DataType getMatrixDataType() const
Definition: LUSolver.h:135
Solves Mx = b for x, where M is a SparseMatrix.
Definition: LUSolver.h:37
Base class for classes that can communicate their status during execution.
Definition: Communicator.h:56
void solve(Matrix< double > &b)
DataType
Definition: LUSolver.h:42