monolish  0.14.2
MONOlithic LIner equation Solvers for Highly-parallel architecture
crs_constructor.cpp
Go to the documentation of this file.
1 #include "../../include/common/monolish_logger.hpp"
2 #include "../../include/common/monolish_matrix.hpp"
3 #include "../../include/common/monolish_vector.hpp"
4 #include "../internal/monolish_internal.hpp"
5 
6 namespace monolish {
7 namespace matrix {
8 
9 // constructor ///
10 template <typename T>
11 CRS<T>::CRS(const size_t M, const size_t N, const size_t NNZ) {
12  Logger &logger = Logger::get_instance();
13  logger.util_in(monolish_func);
14  rowN = M;
15  colN = N;
16  nnz = NNZ;
17  gpu_status = false;
18  row_ptr.resize(M + 1);
19  col_ind.resize(nnz);
20  val.resize(nnz);
21  logger.util_out();
22 }
23 template CRS<double>::CRS(const size_t M, const size_t N, const size_t NNZ);
24 template CRS<float>::CRS(const size_t M, const size_t N, const size_t NNZ);
25 
26 template <typename T>
27 CRS<T>::CRS(const size_t M, const size_t N, const size_t NNZ, const int *rowptr,
28  const int *colind, const T *value) {
29  Logger &logger = Logger::get_instance();
30  logger.util_in(monolish_func);
31  rowN = M;
32  colN = N;
33  nnz = NNZ;
34  gpu_status = false;
35  row_ptr.resize(M + 1);
36  col_ind.resize(nnz);
37  val.resize(nnz);
38  std::copy(rowptr, rowptr + (M + 1), row_ptr.begin());
39  std::copy(colind, colind + nnz, col_ind.begin());
40  std::copy(value, value + nnz, val.begin());
41  compute_hash();
42  logger.util_out();
43 }
44 template CRS<double>::CRS(const size_t M, const size_t N, const size_t NNZ,
45  const int *rowptr, const int *colind,
46  const double *value);
47 template CRS<float>::CRS(const size_t M, const size_t N, const size_t NNZ,
48  const int *rowptr, const int *colind,
49  const float *value);
50 
51 template <typename T>
52 CRS<T>::CRS(const size_t M, const size_t N, const size_t NNZ, const int *rowptr,
53  const int *colind, const T *value, const size_t origin) {
54  Logger &logger = Logger::get_instance();
55  logger.util_in(monolish_func);
56  rowN = M;
57  colN = N;
58  nnz = NNZ;
59  gpu_status = false;
60  row_ptr.resize(M + 1);
61  col_ind.resize(nnz);
62  val.resize(nnz);
63  std::copy(rowptr, rowptr + (M + 1), row_ptr.begin());
64  std::copy(colind, colind + nnz, col_ind.begin());
65  std::copy(value, value + nnz, val.begin());
66 
67 #pragma omp parallel for
68  for (size_t i = 0; i < nnz; i++) {
69  col_ind[i] -= origin;
70  }
71 
72  compute_hash();
73  logger.util_out();
74 }
75 template CRS<double>::CRS(const size_t M, const size_t N, const size_t NNZ,
76  const int *rowptr, const int *colind,
77  const double *value, const size_t origin);
78 template CRS<float>::CRS(const size_t M, const size_t N, const size_t NNZ,
79  const int *rowptr, const int *colind,
80  const float *value, const size_t origin);
81 
82 template <typename T>
83 CRS<T>::CRS(const size_t M, const size_t N, const std::vector<int> &rowptr,
84  const std::vector<int> &colind, const std::vector<T> &value) {
85  Logger &logger = Logger::get_instance();
86  logger.util_in(monolish_func);
87  rowN = M;
88  colN = N;
89  nnz = value.size();
90  gpu_status = false;
91  row_ptr.resize(M + 1);
92  col_ind.resize(nnz);
93  val.resize(nnz);
94  std::copy(rowptr.data(), rowptr.data() + (M + 1), row_ptr.begin());
95  std::copy(colind.data(), colind.data() + nnz, col_ind.begin());
96  std::copy(value.data(), value.data() + nnz, val.begin());
97  compute_hash();
98  logger.util_out();
99 }
100 template CRS<double>::CRS(const size_t M, const size_t N,
101  const std::vector<int> &rowptr,
102  const std::vector<int> &colind,
103  const std::vector<double> &value);
104 template CRS<float>::CRS(const size_t M, const size_t N,
105  const std::vector<int> &rowptr,
106  const std::vector<int> &colind,
107  const std::vector<float> &value);
108 
109 template <typename T>
110 CRS<T>::CRS(const size_t M, const size_t N, const std::vector<int> &rowptr,
111  const std::vector<int> &colind, const vector<T> &value) {
112  Logger &logger = Logger::get_instance();
113  logger.util_in(monolish_func);
114  rowN = M;
115  colN = N;
116  nnz = value.size();
117  gpu_status = false;
118  row_ptr.resize(M + 1);
119  col_ind.resize(nnz);
120  val.resize(nnz);
121  std::copy(rowptr.data(), rowptr.data() + (M + 1), row_ptr.begin());
122  std::copy(colind.data(), colind.data() + nnz, col_ind.begin());
123  std::copy(value.data(), value.data() + nnz, val.begin());
124 
125  if (value.get_device_mem_stat() == true) {
126 #if MONOLISH_USE_NVIDIA_GPU
127  send();
128  const T *data = value.data();
129  T *vald = val.data();
130 #pragma omp target teams distribute parallel for
131  for (size_t i = 0; i < get_nnz(); i++) {
132  vald[i] = data[i];
133  }
134 #else
135  throw std::runtime_error(
136  "error USE_GPU is false, but get_device_mem_stat() == true");
137 #endif
138  }
139 
140  compute_hash();
141  logger.util_out();
142 }
143 
144 template CRS<double>::CRS(const size_t M, const size_t N,
145  const std::vector<int> &rowptr,
146  const std::vector<int> &colind,
147  const vector<double> &value);
148 template CRS<float>::CRS(const size_t M, const size_t N,
149  const std::vector<int> &rowptr,
150  const std::vector<int> &colind,
151  const vector<float> &value);
152 
153 // copy constructor
154 template <typename T> CRS<T>::CRS(const CRS<T> &mat) {
155  Logger &logger = Logger::get_instance();
156  logger.util_in(monolish_func);
157 
158  val.resize(mat.get_nnz());
159  col_ind.resize(mat.get_nnz());
160  row_ptr.resize(mat.get_row() + 1);
161 
162  rowN = mat.get_row();
163  colN = mat.get_col();
164  nnz = mat.get_nnz();
165 
166 #if MONOLISH_USE_NVIDIA_GPU
167  if (mat.get_device_mem_stat()) {
168  send();
169  internal::vcopy(mat.row_ptr.size(), mat.row_ptr.data(), row_ptr.data(),
170  true);
171  internal::vcopy(mat.col_ind.size(), mat.col_ind.data(), col_ind.data(),
172  true);
173  internal::vcopy(mat.val.size(), mat.val.data(), val.data(), true);
174  }
175 #endif
176 
177  internal::vcopy(mat.row_ptr.size(), mat.row_ptr.data(), row_ptr.data(),
178  false);
179  internal::vcopy(mat.col_ind.size(), mat.col_ind.data(), col_ind.data(),
180  false);
181  internal::vcopy(mat.val.size(), mat.val.data(), val.data(), false);
182 
183  compute_hash();
184  logger.util_out();
185 }
186 template CRS<double>::CRS(const CRS<double> &mat);
187 template CRS<float>::CRS(const CRS<float> &mat);
188 } // namespace matrix
189 } // namespace monolish
monolish_func
#define monolish_func
Definition: monolish_logger.hpp:9
monolish::Logger
logger class (singleton, for developper class)
Definition: monolish_logger.hpp:19
monolish::Logger::util_out
void util_out()
Definition: logger_utils.cpp:123
monolish::matrix::CRS::CRS
CRS()
Definition: monolish_crs.hpp:82
monolish::Logger::util_in
void util_in(const std::string func_name)
Definition: logger_utils.cpp:113
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:10
monolish::util::send
auto send(T &x)
send data to GPU
Definition: monolish_common.hpp:642
monolish::Logger::get_instance
static Logger & get_instance()
Definition: monolish_logger.hpp:42