monolish  0.14.0
MONOlithic LIner equation Solvers for Highly-parallel architecture
lu.cpp
Go to the documentation of this file.
1 #include "../../../include/monolish_blas.hpp"
2 #include "../../../include/monolish_equation.hpp"
3 #include "../../internal/lapack/monolish_lapack.hpp"
4 #include "../../internal/monolish_internal.hpp"
5 
6 namespace monolish {
7 
8 template <>
10  vector<double> &x,
11  vector<double> &b) {
12  Logger &logger = Logger::get_instance();
13  logger.func_in(monolish_func);
14 
15  (void)(&A);
16  (void)(&x);
17  (void)(&b);
18  logger.func_out();
19  throw std::runtime_error("error sparse Lu is not impl.");
20 }
21 
22 template <>
24  vector<float> &x,
25  vector<float> &b) {
26  Logger &logger = Logger::get_instance();
27  logger.func_in(monolish_func);
28 
29  (void)(&A);
30  (void)(&x);
31  (void)(&b);
32  logger.func_out();
33  throw std::runtime_error("error sparse Lu is not impl.");
34 }
35 
36 template <>
38  vector<double> &XB) {
39  Logger &logger = Logger::get_instance();
40  logger.func_in(monolish_func);
41 
42  int ret = MONOLISH_SOLVER_SUCCESS;
43 
44  if (lib == 1) {
45  std::vector<int> ipiv(std::min(A.get_row(), A.get_col()));
46 
47  if (internal::lapack::getrf(A, ipiv) != 0) {
49  };
50 
51  if (internal::lapack::getrs(A, XB, ipiv) != 0) {
53  }
54 
55  } else {
56  logger.func_out();
57  throw std::runtime_error("error solver.lib is not 1");
58  }
59 
60  logger.func_out();
61  return ret;
62 }
63 template <>
65  vector<float> &XB) {
66  Logger &logger = Logger::get_instance();
67  logger.func_in(monolish_func);
68 
69  int ret = MONOLISH_SOLVER_SUCCESS;
70 
71  if (lib == 1) {
72  std::vector<int> ipiv(std::min(A.get_row(), A.get_col()));
73 
74  if (internal::lapack::getrf(A, ipiv) != 0) {
76  };
77 
78  if (internal::lapack::getrs(A, XB, ipiv) != 0) {
80  }
81  } else {
82  logger.func_out();
83  throw std::runtime_error("error solver.lib is not 1");
84  }
85 
86  logger.func_out();
87  return ret;
88 }
89 
90 template <>
92  vector<double> &x,
93  vector<double> &b) {
94  Logger &logger = Logger::get_instance();
95  logger.func_in(monolish_func);
96 
97  int ret = MONOLISH_SOLVER_SUCCESS;
98 
99  if (lib == 1) {
100  std::vector<int> ipiv(std::min(A.get_row(), A.get_col()));
101  monolish::blas::copy(b, x);
102 
103  if (internal::lapack::getrf(A, ipiv) != 0) {
105  };
106 
107  if (internal::lapack::getrs(A, x, ipiv) != 0) {
109  }
110 
111  } else {
112  logger.func_out();
113  throw std::runtime_error("error solver.lib is not 1");
114  }
115 
116  logger.func_out();
117  return ret;
118 }
119 
120 template <>
122  vector<float> &x,
123  vector<float> &b) {
124  Logger &logger = Logger::get_instance();
125  logger.func_in(monolish_func);
126 
127  int ret = MONOLISH_SOLVER_SUCCESS;
128 
129  if (lib == 1) {
130  std::vector<int> ipiv(std::min(A.get_row(), A.get_col()));
131  monolish::blas::copy(b, x);
132 
133  if (internal::lapack::getrf(A, ipiv) != 0) {
135  };
136 
137  if (internal::lapack::getrs(A, x, ipiv) != 0) {
139  }
140  } else {
141  logger.func_out();
142  throw std::runtime_error("error solver.lib is not 1");
143  }
144 
145  logger.func_out();
146  return ret;
147 }
148 
149 } // namespace monolish
monolish_func
#define monolish_func
Definition: monolish_logger.hpp:9
monolish::vml::min
void min(const matrix::CRS< double > &A, const matrix::CRS< double > &B, matrix::CRS< double > &C)
Create a new CRS matrix with smallest elements of two matrices (C[0:nnz] = min(A[0:nnz],...
Definition: matrix_vml.cpp:390
monolish::Logger
logger class (singleton, for developper class)
Definition: monolish_logger.hpp:19
monolish::Logger::func_out
void func_out()
Definition: logger_utils.cpp:80
monolish::matrix::Dense::get_row
size_t get_row() const
get # of row
Definition: monolish_dense.hpp:199
monolish::matrix::Dense
Dense format Matrix.
Definition: monolish_coo.hpp:28
MONOLISH_SOLVER_SUCCESS
#define MONOLISH_SOLVER_SUCCESS
Definition: monolish_common.hpp:10
MONOLISH_SOLVER_BREAKDOWN
#define MONOLISH_SOLVER_BREAKDOWN
Definition: monolish_common.hpp:13
monolish::equation::LU
LU solver class (Dense, CPU only now)
Definition: monolish_equation.hpp:201
monolish::blas::copy
void copy(const matrix::Dense< double > &A, matrix::Dense< double > &C)
Dense matrix copy (C=A)
Definition: dense_copy.cpp:25
monolish
Definition: monolish_matrix_blas.hpp:9
monolish::matrix::Dense::get_col
size_t get_col() const
get # of col
Definition: monolish_dense.hpp:208
monolish::vector
vector class
Definition: monolish_coo.hpp:25
monolish::Logger::get_instance
static Logger & get_instance()
Definition: monolish_logger.hpp:42
monolish::matrix::CRS
Compressed Row Storage (CRS) format Matrix.
Definition: monolish_coo.hpp:29
monolish::Logger::func_in
void func_in(const std::string func_name)
Definition: logger_utils.cpp:69