ViennaCL - The Vienna Computing Library  1.7.1
Free open-source GPU-accelerated linear algebra and solver library.
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
misc_operations.hpp
Go to the documentation of this file.
1 #ifndef VIENNACL_LINALG_MISC_OPERATIONS_HPP_
2 #define VIENNACL_LINALG_MISC_OPERATIONS_HPP_
3 
4 /* =========================================================================
5  Copyright (c) 2010-2016, 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/matrix.hpp"
29 #include "viennacl/tools/tools.hpp"
31 
32 #ifdef VIENNACL_WITH_OPENCL
34 #endif
35 
36 #ifdef VIENNACL_WITH_CUDA
38 #endif
39 
40 namespace viennacl
41 {
42  namespace linalg
43  {
44 
45  namespace detail
46  {
47 
48  template<typename ScalarType>
50  viennacl::backend::mem_handle const & row_index_array,
51  viennacl::backend::mem_handle const & row_buffer,
52  viennacl::backend::mem_handle const & col_buffer,
53  viennacl::backend::mem_handle const & element_buffer,
54  vcl_size_t num_rows
55  )
56  {
57  assert( viennacl::traits::handle(vec).get_active_handle_id() == row_index_array.get_active_handle_id() && bool("Incompatible memory domains"));
58  assert( viennacl::traits::handle(vec).get_active_handle_id() == row_buffer.get_active_handle_id() && bool("Incompatible memory domains"));
59  assert( viennacl::traits::handle(vec).get_active_handle_id() == col_buffer.get_active_handle_id() && bool("Incompatible memory domains"));
60  assert( viennacl::traits::handle(vec).get_active_handle_id() == element_buffer.get_active_handle_id() && bool("Incompatible memory domains"));
61 
63  {
65  viennacl::linalg::host_based::detail::level_scheduling_substitute(vec, row_index_array, row_buffer, col_buffer, element_buffer, num_rows);
66  break;
67 #ifdef VIENNACL_WITH_OPENCL
69  viennacl::linalg::opencl::detail::level_scheduling_substitute(vec, row_index_array, row_buffer, col_buffer, element_buffer, num_rows);
70  break;
71 #endif
72 #ifdef VIENNACL_WITH_CUDA
74  viennacl::linalg::cuda::detail::level_scheduling_substitute(vec, row_index_array, row_buffer, col_buffer, element_buffer, num_rows);
75  break;
76 #endif
78  throw memory_exception("not initialised!");
79  default:
80  throw memory_exception("not implemented");
81  }
82  }
83 
84 
85 
86 
87  } //namespace detail
88 
89 
90  } //namespace linalg
91 } //namespace viennacl
92 
93 
94 #endif
Exception class in case of memory errors.
Definition: forwards.h:572
Implementation of the dense matrix class.
Various little tools used here and there in ViennaCL.
This file provides the forward declarations for the main types used within ViennaCL.
void level_scheduling_substitute(viennacl::vector< NumericT > &vec, std::list< viennacl::backend::mem_handle > const &row_index_arrays, std::list< viennacl::backend::mem_handle > const &row_buffers, std::list< viennacl::backend::mem_handle > const &col_buffers, std::list< viennacl::backend::mem_handle > const &element_buffers, std::list< vcl_size_t > const &row_elimination_num_list)
Definition: common.hpp:224
Implementations of miscellaneous operations on the CPU using a single thread or OpenMP.
Implementations of operations using compressed_matrix and OpenCL.
std::size_t vcl_size_t
Definition: forwards.h:75
void level_scheduling_substitute(vector< NumericT > &x, 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)
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
viennacl::backend::mem_handle & handle(T &obj)
Returns the generic memory handle of an object. Non-const version.
Definition: handle.hpp:41
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)
Implementations of miscellaneous operations using CUDA.
Implementation of the ViennaCL scalar class.
memory_types get_active_handle_id() const
Returns an ID for the currently active memory buffer. Other memory buffers might contain old or no da...
Definition: mem_handle.hpp:118