monolish  0.14.2
MONOlithic LIner equation Solvers for Highly-parallel architecture
nrm2.hpp
Go to the documentation of this file.
1 #pragma once
2 
3 namespace monolish {
4 
5 namespace {
6 template <typename F1> double Dnrm2_core(const F1 &x) {
7  Logger &logger = Logger::get_instance();
8  logger.func_in(monolish_func);
9 
10  double ans = 0;
11  const double *xd = x.data();
12  size_t size = x.size();
13  const size_t xoffset = x.get_offset();
14 
15  if (x.get_device_mem_stat() == true) {
16 #if MONOLISH_USE_NVIDIA_GPU
17  cublasHandle_t h;
18  internal::check_CUDA(cublasCreate(&h));
19 #pragma omp target data use_device_ptr(xd)
20  {
21  internal::check_CUDA(
22  cublasDdot(h, size, xd + xoffset, 1, xd + xoffset, 1, &ans));
23  }
24  cublasDestroy(h);
25 #else
26  throw std::runtime_error(
27  "error USE_GPU is false, but get_device_mem_stat() == true");
28 #endif
29  } else {
30  ans = cblas_ddot(size, xd + xoffset, 1, xd + xoffset, 1);
31  }
32 
33 #if MONOLISH_USE_MPI
34  mpi::comm &comm = mpi::comm::get_instance();
35  ans = comm.Allreduce(ans);
36 #endif
37 
38  logger.func_out();
39  return sqrt(ans);
40 }
41 
42 template <typename F1> float Snrm2_core(const F1 &x) {
43  Logger &logger = Logger::get_instance();
44  logger.func_in(monolish_func);
45 
46  float ans = 0;
47  const float *xd = x.data();
48  size_t size = x.size();
49  const size_t xoffset = x.get_offset();
50 
51  if (x.get_device_mem_stat() == true) {
52 #if MONOLISH_USE_NVIDIA_GPU
53  cublasHandle_t h;
54  internal::check_CUDA(cublasCreate(&h));
55 #pragma omp target data use_device_ptr(xd)
56  {
57  internal::check_CUDA(
58  cublasSdot(h, size, xd + xoffset, 1, xd + xoffset, 1, &ans));
59  }
60  cublasDestroy(h);
61 #else
62  throw std::runtime_error(
63  "error USE_GPU is false, but get_device_mem_stat() == true");
64 #endif
65  } else {
66  ans = cblas_sdot(size, xd + xoffset, 1, xd + xoffset, 1);
67  }
68 
69 #if MONOLISH_USE_MPI
70  mpi::comm &comm = mpi::comm::get_instance();
71  ans = comm.Allreduce(ans);
72 #endif
73 
74  logger.func_out();
75  return sqrt(ans);
76 }
77 } // namespace
78 
79 } // namespace monolish
monolish_func
#define monolish_func
Definition: monolish_logger.hpp:9
monolish::vml::sqrt
void sqrt(const matrix::CRS< double > &A, matrix::CRS< double > &C)
sqrt to CRS matrix elements (C[0:nnz] = sqrt(A[0:nnz]))
Definition: matrix_vml.cpp:309
monolish
Definition: monolish_matrix_blas.hpp:10
monolish::mpi::comm::get_instance
static comm & get_instance()
Definition: monolish_mpi_core.hpp:33
monolish::Logger::get_instance
static Logger & get_instance()
Definition: monolish_logger.hpp:42