monolish  0.14.0
MONOlithic LIner equation Solvers for Highly-parallel architecture
dense_getvec.cpp
Go to the documentation of this file.
1 #include "../../../include/monolish_blas.hpp"
2 #include "../../internal/monolish_internal.hpp"
3 
4 namespace monolish {
5 namespace matrix {
6 
7 // diag
8 template <typename T> void Dense<T>::diag(vector<T> &vec) const {
9  Logger &logger = Logger::get_instance();
10  logger.func_in(monolish_func);
11 
12  T *vecd = vec.data();
13 
14  const T *vald = val.data();
15  const size_t N = get_col();
16  const size_t Len = std::min(get_row(), get_col());
17 
18  assert(Len == vec.size());
19  assert(get_device_mem_stat() == vec.get_device_mem_stat());
20 
21  if (gpu_status == true) {
22 #if MONOLISH_USE_GPU // gpu
23 #pragma omp target teams distribute parallel for
24  for (size_t i = 0; i < Len; i++) {
25  vecd[i] = vald[N * i + i];
26  }
27 #else
28  throw std::runtime_error("error USE_GPU is false, but gpu_status == true");
29 #endif
30  } else {
31 #pragma omp parallel for
32  for (size_t i = 0; i < Len; i++) {
33  vecd[i] = vald[N * i + i];
34  }
35  }
36 
37  logger.func_out();
38 }
39 template void monolish::matrix::Dense<double>::diag(vector<double> &vec) const;
40 template void monolish::matrix::Dense<float>::diag(vector<float> &vec) const;
41 
42 template <typename T> void Dense<T>::diag(view1D<vector<T>, T> &vec) const {
43  Logger &logger = Logger::get_instance();
44  logger.func_in(monolish_func);
45 
46  T *vecd = vec.data();
47 
48  const T *vald = val.data();
49  const size_t N = get_col();
50  const size_t Len = std::min(get_row(), get_col());
51 
52  assert(Len == vec.size());
53  assert(get_device_mem_stat() == vec.get_device_mem_stat());
54 
55  if (gpu_status == true) {
56 #if MONOLISH_USE_GPU // gpu
57 #pragma omp target teams distribute parallel for
58  for (size_t i = 0; i < Len; i++) {
59  vecd[i] = vald[N * i + i];
60  }
61 #else
62  throw std::runtime_error("error USE_GPU is false, but gpu_status == true");
63 #endif
64  } else {
65 #pragma omp parallel for
66  for (size_t i = 0; i < Len; i++) {
67  vecd[i] = vald[N * i + i];
68  }
69  }
70 
71  logger.func_out();
72 }
74  view1D<vector<double>, double> &vec) const;
75 template void
76 monolish::matrix::Dense<float>::diag(view1D<vector<float>, float> &vec) const;
77 
78 template <typename T>
79 void Dense<T>::diag(view1D<matrix::Dense<T>, T> &vec) const {
80  Logger &logger = Logger::get_instance();
81  logger.func_in(monolish_func);
82 
83  T *vecd = vec.data();
84 
85  const T *vald = val.data();
86  const size_t N = get_col();
87  const size_t Len = std::min(get_row(), get_col());
88 
89  assert(Len == vec.size());
90  assert(get_device_mem_stat() == vec.get_device_mem_stat());
91 
92  if (gpu_status == true) {
93 #if MONOLISH_USE_GPU // gpu
94 #pragma omp target teams distribute parallel for
95  for (size_t i = 0; i < Len; i++) {
96  vecd[i] = vald[N * i + i];
97  }
98 #else
99  throw std::runtime_error("error USE_GPU is false, but gpu_status == true");
100 #endif
101  } else {
102 #pragma omp parallel for
103  for (size_t i = 0; i < Len; i++) {
104  vecd[i] = vald[N * i + i];
105  }
106  }
107 
108  logger.func_out();
109 }
111  view1D<matrix::Dense<double>, double> &vec) const;
113  view1D<matrix::Dense<float>, float> &vec) const;
114 
115 // row
116 template <typename T> void Dense<T>::row(const size_t r, vector<T> &vec) const {
117  Logger &logger = Logger::get_instance();
118  logger.func_in(monolish_func);
119 
120  T *vecd = vec.data();
121 
122  const T *vald = val.data();
123  const size_t N = get_col();
124 
125  assert(N == vec.size());
126  assert(get_device_mem_stat() == vec.get_device_mem_stat());
127 
128  if (gpu_status == true) {
129 #if MONOLISH_USE_GPU // gpu
130 #pragma omp target teams distribute parallel for
131  for (size_t i = 0; i < N; i++) {
132  vecd[i] = vald[r * N + i];
133  }
134 #else
135  throw std::runtime_error("error USE_GPU is false, but gpu_status == true");
136 #endif
137  } else {
138 #pragma omp parallel for
139  for (size_t i = 0; i < N; i++) {
140  vecd[i] = vald[r * N + i];
141  }
142  }
143 
144  logger.func_out();
145 }
146 template void monolish::matrix::Dense<double>::row(const size_t r,
147  vector<double> &vec) const;
148 template void monolish::matrix::Dense<float>::row(const size_t r,
149  vector<float> &vec) const;
150 
151 template <typename T>
152 void Dense<T>::row(const size_t r, view1D<vector<T>, T> &vec) const {
153  Logger &logger = Logger::get_instance();
154  logger.func_in(monolish_func);
155 
156  T *vecd = vec.data();
157 
158  const T *vald = val.data();
159  const size_t N = get_col();
160 
161  assert(N == vec.size());
162  assert(get_device_mem_stat() == vec.get_device_mem_stat());
163 
164  if (gpu_status == true) {
165 #if MONOLISH_USE_GPU // gpu
166 #pragma omp target teams distribute parallel for
167  for (size_t i = 0; i < N; i++) {
168  vecd[i] = vald[r * N + i];
169  }
170 #else
171  throw std::runtime_error("error USE_GPU is false, but gpu_status == true");
172 #endif
173  } else {
174 #pragma omp parallel for
175  for (size_t i = 0; i < N; i++) {
176  vecd[i] = vald[r * N + i];
177  }
178  }
179 
180  logger.func_out();
181 }
182 template void
184  view1D<vector<double>, double> &vec) const;
185 template void
187  view1D<vector<float>, float> &vec) const;
188 
189 template <typename T>
190 void Dense<T>::row(const size_t r, view1D<matrix::Dense<T>, T> &vec) const {
191  Logger &logger = Logger::get_instance();
192  logger.func_in(monolish_func);
193 
194  T *vecd = vec.data();
195 
196  const T *vald = val.data();
197  const size_t N = get_col();
198 
199  assert(N == vec.size());
200  assert(get_device_mem_stat() == vec.get_device_mem_stat());
201 
202  if (gpu_status == true) {
203 #if MONOLISH_USE_GPU // gpu
204 #pragma omp target teams distribute parallel for
205  for (size_t i = 0; i < N; i++) {
206  vecd[i] = vald[r * N + i];
207  }
208 #else
209  throw std::runtime_error("error USE_GPU is false, but gpu_status == true");
210 #endif
211  } else {
212 #pragma omp parallel for
213  for (size_t i = 0; i < N; i++) {
214  vecd[i] = vald[r * N + i];
215  }
216  }
217 
218  logger.func_out();
219 }
221  const size_t r, view1D<matrix::Dense<double>, double> &vec) const;
223  const size_t r, view1D<matrix::Dense<float>, float> &vec) const;
224 
225 // col
226 template <typename T> void Dense<T>::col(const size_t c, vector<T> &vec) const {
227  Logger &logger = Logger::get_instance();
228  logger.func_in(monolish_func);
229 
230  T *vecd = vec.data();
231 
232  const T *vald = val.data();
233  const size_t M = get_row();
234  const size_t N = get_col();
235 
236  assert(M == vec.size());
237  assert(get_device_mem_stat() == vec.get_device_mem_stat());
238 
239  if (gpu_status == true) {
240 #if MONOLISH_USE_GPU // gpu
241 #pragma omp target teams distribute parallel for
242  for (size_t i = 0; i < M; i++) {
243  vecd[i] = vald[i * N + c];
244  }
245 #else
246  throw std::runtime_error("error USE_GPU is false, but gpu_status == true");
247 #endif
248  } else {
249 #pragma omp parallel for
250  for (size_t i = 0; i < M; i++) {
251  vecd[i] = vald[i * N + c];
252  }
253  }
254 
255  logger.func_out();
256 }
257 template void monolish::matrix::Dense<double>::col(const size_t c,
258  vector<double> &vec) const;
259 template void monolish::matrix::Dense<float>::col(const size_t c,
260  vector<float> &vec) const;
261 
262 template <typename T>
263 void Dense<T>::col(const size_t c, view1D<vector<T>, T> &vec) const {
264  Logger &logger = Logger::get_instance();
265  logger.func_in(monolish_func);
266 
267  T *vecd = vec.data();
268 
269  const T *vald = val.data();
270  const size_t M = get_row();
271  const size_t N = get_col();
272 
273  assert(M == vec.size());
274  assert(get_device_mem_stat() == vec.get_device_mem_stat());
275 
276  if (gpu_status == true) {
277 #if MONOLISH_USE_GPU // gpu
278 #pragma omp target teams distribute parallel for
279  for (size_t i = 0; i < M; i++) {
280  vecd[i] = vald[i * N + c];
281  }
282 #else
283  throw std::runtime_error("error USE_GPU is false, but gpu_status == true");
284 #endif
285  } else {
286 #pragma omp parallel for
287  for (size_t i = 0; i < M; i++) {
288  vecd[i] = vald[i * N + c];
289  }
290  }
291 
292  logger.func_out();
293 }
294 template void
296  view1D<vector<double>, double> &vec) const;
297 template void
299  view1D<vector<float>, float> &vec) const;
300 
301 template <typename T>
302 void Dense<T>::col(const size_t c, view1D<matrix::Dense<T>, T> &vec) const {
303  Logger &logger = Logger::get_instance();
304  logger.func_in(monolish_func);
305 
306  T *vecd = vec.data();
307 
308  const T *vald = val.data();
309  const size_t M = get_row();
310  const size_t N = get_col();
311 
312  assert(M == vec.size());
313  assert(get_device_mem_stat() == vec.get_device_mem_stat());
314 
315  if (gpu_status == true) {
316 #if MONOLISH_USE_GPU // gpu
317 #pragma omp target teams distribute parallel for
318  for (size_t i = 0; i < M; i++) {
319  vecd[i] = vald[i * N + c];
320  }
321 #else
322  throw std::runtime_error("error USE_GPU is false, but gpu_status == true");
323 #endif
324  } else {
325 #pragma omp parallel for
326  for (size_t i = 0; i < M; i++) {
327  vecd[i] = vald[i * N + c];
328  }
329  }
330 
331  logger.func_out();
332 }
334  const size_t c, view1D<matrix::Dense<double>, double> &vec) const;
336  const size_t c, view1D<matrix::Dense<float>, float> &vec) const;
337 } // namespace matrix
338 } // 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::matrix::Dense::col
void col(const size_t c, vector< Float > &vec) const
get column vector
monolish
Definition: monolish_matrix_blas.hpp:9
monolish::Logger::get_instance
static Logger & get_instance()
Definition: monolish_logger.hpp:42
monolish::matrix::Dense::row
void row(const size_t r, vector< Float > &vec) const
get row vector
monolish::matrix::Dense::diag
void diag(vector< Float > &vec) const
get diag. vector