ViennaCL - The Vienna Computing Library  1.7.0
Free open-source GPU-accelerated linear algebra and solver library.
misc_operations.hpp
Go to the documentation of this file.
1 #ifndef VIENNACL_LINALG_CUDA_MISC_OPERATIONS_HPP_
2 #define VIENNACL_LINALG_CUDA_MISC_OPERATIONS_HPP_
3 
4 /* =========================================================================
5  Copyright (c) 2010-2015, Institute for Microelectronics,
6  Institute for Analysis and Scientific Computing,
7  TU Wien.
8  Portions of this software are copyright by UChicago Argonne, LLC.
9 
10  -----------------
11  ViennaCL - The Vienna Computing Library
12  -----------------
13 
14  Project Head: Karl Rupp rupp@iue.tuwien.ac.at
15 
16  (A list of authors and contributors can be found in the manual)
17 
18  License: MIT (X11), see file LICENSE in the base directory
19 ============================================================================= */
20 
25 #include "viennacl/forwards.h"
26 #include "viennacl/scalar.hpp"
27 #include "viennacl/vector.hpp"
28 #include "viennacl/tools/tools.hpp"
30 
31 
32 namespace viennacl
33 {
34 namespace linalg
35 {
36 namespace cuda
37 {
38 namespace detail
39 {
40 
41 template<typename NumericT>
43  const unsigned int * row_index_array,
44  const unsigned int * row_indices,
45  const unsigned int * column_indices,
46  const NumericT * elements,
47  NumericT * vec,
48  unsigned int size)
49 {
50  for (unsigned int row = blockDim.x * blockIdx.x + threadIdx.x;
51  row < size;
52  row += gridDim.x * blockDim.x)
53  {
54  unsigned int eq_row = row_index_array[row];
55  NumericT vec_entry = vec[eq_row];
56  unsigned int row_end = row_indices[row+1];
57 
58  for (unsigned int j = row_indices[row]; j < row_end; ++j)
59  vec_entry -= vec[column_indices[j]] * elements[j];
60 
61  vec[eq_row] = vec_entry;
62  }
63 }
64 
65 
66 
67 template<typename NumericT>
69  viennacl::backend::mem_handle const & row_index_array,
70  viennacl::backend::mem_handle const & row_buffer,
71  viennacl::backend::mem_handle const & col_buffer,
72  viennacl::backend::mem_handle const & element_buffer,
73  vcl_size_t num_rows
74  )
75 {
76  level_scheduling_substitute_kernel<<<128, 128>>>(viennacl::cuda_arg<unsigned int>(row_index_array),
77  viennacl::cuda_arg<unsigned int>(row_buffer),
78  viennacl::cuda_arg<unsigned int>(col_buffer),
79  viennacl::cuda_arg<NumericT>(element_buffer),
80  viennacl::cuda_arg(vec),
81  static_cast<unsigned int>(num_rows)
82  );
83 }
84 
85 } //namespace detail
86 } //namespace cuda
87 } //namespace linalg
88 } //namespace viennacl
89 
90 
91 #endif
Various little tools used here and there in ViennaCL.
This file provides the forward declarations for the main types used within ViennaCL.
float NumericT
Definition: bisect.cpp:40
Main namespace in ViennaCL. Holds all the basic types such as vector, matrix, etc. and defines operations upon them.
Definition: cpu_ram.hpp:34
vcl_size_t size(VectorType const &vec)
Generic routine for obtaining the size of a vector (ViennaCL, uBLAS, etc.)
Definition: size.hpp:235
Definition: blas3.hpp:36
std::size_t vcl_size_t
Definition: forwards.h:75
vector_expression< const matrix_base< NumericT, F >, const unsigned int, op_row > row(const matrix_base< NumericT, F > &A, unsigned int i)
Definition: matrix.hpp:900
Common routines for CUDA execution.
The vector type with operator-overloads and proxy classes is defined here. Linear algebra operations ...
void level_scheduling_substitute(vector< NumericT > &vec, viennacl::backend::mem_handle const &row_index_array, viennacl::backend::mem_handle const &row_buffer, viennacl::backend::mem_handle const &col_buffer, viennacl::backend::mem_handle const &element_buffer, vcl_size_t num_rows)
Main abstraction class for multiple memory domains. Represents a buffer in either main RAM...
Definition: mem_handle.hpp:89
__global__ void level_scheduling_substitute_kernel(const unsigned int *row_index_array, const unsigned int *row_indices, const unsigned int *column_indices, const NumericT *elements, NumericT *vec, unsigned int size)
NumericT * cuda_arg(scalar< NumericT > &obj)
Convenience helper function for extracting the CUDA handle from a ViennaCL scalar. Non-const version.
Definition: common.hpp:39
Implementation of the ViennaCL scalar class.