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
ichol.hpp
Go to the documentation of this file.
1 #ifndef VIENNACL_LINALG_ICHOL_HPP_
2 #define VIENNACL_LINALG_ICHOL_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 <vector>
26 #include <cmath>
27 #include <iostream>
28 #include "viennacl/forwards.h"
29 #include "viennacl/tools/tools.hpp"
31 
33 
34 #include <map>
35 
36 namespace viennacl
37 {
38 namespace linalg
39 {
40 
43 class ichol0_tag {};
44 
45 
54 template<typename NumericT>
56 {
57  assert( (viennacl::traits::context(A).memory_type() == viennacl::MAIN_MEMORY) && bool("System matrix must reside in main memory for ICHOL0") );
58 
59  NumericT * elements = viennacl::linalg::host_based::detail::extract_raw_pointer<NumericT>(A.handle());
60  unsigned int const * row_buffer = viennacl::linalg::host_based::detail::extract_raw_pointer<unsigned int>(A.handle1());
61  unsigned int const * col_buffer = viennacl::linalg::host_based::detail::extract_raw_pointer<unsigned int>(A.handle2());
62 
63  //std::cout << A.size1() << std::endl;
64  for (vcl_size_t i=0; i<A.size1(); ++i)
65  {
66  unsigned int row_i_begin = row_buffer[i];
67  unsigned int row_i_end = row_buffer[i+1];
68 
69  // get a_ii:
70  NumericT a_ii = 0;
71  for (unsigned int buf_index_aii = row_i_begin; buf_index_aii < row_i_end; ++buf_index_aii)
72  {
73  if (col_buffer[buf_index_aii] == i)
74  {
75  a_ii = std::sqrt(elements[buf_index_aii]);
76  elements[buf_index_aii] = a_ii;
77  break;
78  }
79  }
80 
81  // Now scale column/row i, i.e. A(k, i) /= A(i, i)
82  for (unsigned int buf_index_aii = row_i_begin; buf_index_aii < row_i_end; ++buf_index_aii)
83  {
84  if (col_buffer[buf_index_aii] > i)
85  elements[buf_index_aii] /= a_ii;
86  }
87 
88  // Now compute A(k, j) -= A(k, i) * A(j, i) for all nonzero k, j in column i:
89  for (unsigned int buf_index_j = row_i_begin; buf_index_j < row_i_end; ++buf_index_j)
90  {
91  unsigned int j = col_buffer[buf_index_j];
92  if (j <= i)
93  continue;
94 
95  NumericT a_ji = elements[buf_index_j];
96 
97  for (unsigned int buf_index_k = row_i_begin; buf_index_k < row_i_end; ++buf_index_k)
98  {
99  unsigned int k = col_buffer[buf_index_k];
100  if (k < j)
101  continue;
102 
103  NumericT a_ki = elements[buf_index_k];
104 
105  //Now check whether A(k, j) is in nonzero pattern:
106  unsigned int row_j_begin = row_buffer[j];
107  unsigned int row_j_end = row_buffer[j+1];
108  for (unsigned int buf_index_kj = row_j_begin; buf_index_kj < row_j_end; ++buf_index_kj)
109  {
110  if (col_buffer[buf_index_kj] == k)
111  {
112  elements[buf_index_kj] -= a_ki * a_ji;
113  break;
114  }
115  }
116  }
117  }
118 
119  }
120 
121 }
122 
123 
126 template<typename MatrixT>
128 {
129  typedef typename MatrixT::value_type NumericType;
130 
131 public:
132  ichol0_precond(MatrixT const & mat, ichol0_tag const & tag) : tag_(tag), LLT(mat.size1(), mat.size2(), viennacl::context(viennacl::MAIN_MEMORY))
133  {
134  //initialize preconditioner:
135  //std::cout << "Start CPU precond" << std::endl;
136  init(mat);
137  //std::cout << "End CPU precond" << std::endl;
138  }
139 
140  template<typename VectorT>
141  void apply(VectorT & vec) const
142  {
143  unsigned int const * row_buffer = viennacl::linalg::host_based::detail::extract_raw_pointer<unsigned int>(LLT.handle1());
144  unsigned int const * col_buffer = viennacl::linalg::host_based::detail::extract_raw_pointer<unsigned int>(LLT.handle2());
145  NumericType const * elements = viennacl::linalg::host_based::detail::extract_raw_pointer<NumericType>(LLT.handle());
146 
147  // Note: L is stored in a column-oriented fashion, i.e. transposed w.r.t. the row-oriented layout. Thus, the factorization A = L L^T holds L in the upper triangular part of A.
148  viennacl::linalg::host_based::detail::csr_trans_inplace_solve<NumericType>(row_buffer, col_buffer, elements, vec, LLT.size2(), lower_tag());
149  viennacl::linalg::host_based::detail::csr_inplace_solve<NumericType>(row_buffer, col_buffer, elements, vec, LLT.size2(), upper_tag());
150  }
151 
152 private:
153  void init(MatrixT const & mat)
154  {
156  viennacl::switch_memory_context(LLT, host_ctx);
157 
158  viennacl::copy(mat, LLT);
160  }
161 
162  ichol0_tag const & tag_;
164 };
165 
166 
171 template<typename NumericT, unsigned int AlignmentV>
173 {
175 
176 public:
177  ichol0_precond(MatrixType const & mat, ichol0_tag const & tag) : tag_(tag), LLT(mat.size1(), mat.size2(), viennacl::traits::context(mat))
178  {
179  //initialize preconditioner:
180  //std::cout << "Start GPU precond" << std::endl;
181  init(mat);
182  //std::cout << "End GPU precond" << std::endl;
183  }
184 
185  void apply(vector<NumericT> & vec) const
186  {
187  if (viennacl::traits::context(vec).memory_type() != viennacl::MAIN_MEMORY)
188  {
191 
192  viennacl::switch_memory_context(vec, host_ctx);
195  viennacl::switch_memory_context(vec, old_ctx);
196  }
197  else //apply ILU0 directly:
198  {
199  // Note: L is stored in a column-oriented fashion, i.e. transposed w.r.t. the row-oriented layout. Thus, the factorization A = L L^T holds L in the upper triangular part of A.
202  }
203  }
204 
205 private:
206  void init(MatrixType const & mat)
207  {
209  viennacl::switch_memory_context(LLT, host_ctx);
210  LLT = mat;
211 
213  }
214 
215  ichol0_tag const & tag_;
217 };
218 
219 }
220 }
221 
222 
223 
224 
225 #endif
226 
227 
228 
const vcl_size_t & size2() const
Returns the number of columns.
void inplace_solve(const matrix_base< NumericT > &A, matrix_base< NumericT > &B, SolverTagT)
Direct inplace solver for triangular systems with multiple right hand sides, i.e. A \ B (MATLAB notat...
Incomplete Cholesky preconditioner class with static pattern (ICHOL0), can be supplied to solve()-rou...
Definition: ichol.hpp:127
A tag for incomplete Cholesky factorization with static pattern (ILU0)
Definition: ichol.hpp:43
const vcl_size_t & size1() const
Returns the number of rows.
Various little tools used here and there in ViennaCL.
vcl_size_t size1(MatrixType const &mat)
Generic routine for obtaining the number of rows of a matrix (ViennaCL, uBLAS, etc.)
Definition: size.hpp:163
A tag class representing a lower triangular matrix.
Definition: forwards.h:849
void precondition(viennacl::compressed_matrix< NumericT > &A, ilu0_tag const &)
Implementation of a ILU-preconditioner with static pattern. Optimized version for CSR matrices...
Definition: ilu0.hpp:78
This file provides the forward declarations for the main types used within ViennaCL.
result_of::size_type< MatrixType >::type size2(MatrixType const &mat)
Generic routine for obtaining the number of columns of a matrix (ViennaCL, uBLAS, etc...
Definition: size.hpp:201
const handle_type & handle() const
Returns the OpenCL handle to the matrix entry array.
const handle_type & handle1() const
Returns the OpenCL handle to the row index array.
float NumericT
Definition: bisect.cpp:40
Represents a generic 'context' similar to an OpenCL context, but is backend-agnostic and thus also su...
Definition: context.hpp:39
A tag class representing an upper triangular matrix.
Definition: forwards.h:854
Implementation of the compressed_matrix class.
const handle_type & handle2() const
Returns the OpenCL handle to the column index array.
std::size_t vcl_size_t
Definition: forwards.h:75
ichol0_precond(MatrixT const &mat, ichol0_tag const &tag)
Definition: ichol.hpp:132
Common routines for single-threaded or OpenMP-enabled execution on CPU.
void trans(const matrix_expression< const matrix_base< NumericT, SizeT, DistanceT >, const matrix_base< NumericT, SizeT, DistanceT >, op_trans > &proxy, matrix_base< NumericT > &temp_trans)
void apply(VectorT &vec) const
Definition: ichol.hpp:141
viennacl::context context(T const &t)
Returns an ID for the currently active memory domain of an object.
Definition: context.hpp:40
void copy(std::vector< NumericT > &cpu_vec, circulant_matrix< NumericT, AlignmentV > &gpu_mat)
Copies a circulant matrix from the std::vector to the OpenCL device (either GPU or multi-core CPU) ...
ichol0_precond(MatrixType const &mat, ichol0_tag const &tag)
Definition: ichol.hpp:177
void switch_memory_context(T &obj, viennacl::context new_ctx)
Generic convenience routine for migrating data of an object to a new memory domain.
Definition: memory.hpp:622