ViennaCL - The Vienna Computing Library  1.7.0
Free open-source GPU-accelerated linear algebra and solver library.
ilu_operations.hpp
Go to the documentation of this file.
1 #ifndef VIENNACL_LINALG_ILU_OPERATIONS_HPP_
2 #define VIENNACL_LINALG_ILU_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 PDF manual)
17 
18  License: MIT (X11), see file LICENSE in the base directory
19 ============================================================================= */
20 
25 #include "viennacl/forwards.h"
26 #include "viennacl/range.hpp"
27 #include "viennacl/scalar.hpp"
28 #include "viennacl/tools/tools.hpp"
31 #include "viennacl/traits/size.hpp"
36 
37 #ifdef VIENNACL_WITH_OPENCL
39 #endif
40 
41 #ifdef VIENNACL_WITH_CUDA
43 #endif
44 
45 namespace viennacl
46 {
47 namespace linalg
48 {
49 
55 template<typename NumericT>
58 {
59  switch (viennacl::traits::handle(A).get_active_handle_id())
60  {
63  break;
64 #ifdef VIENNACL_WITH_OPENCL
67  break;
68 #endif
69 #ifdef VIENNACL_WITH_CUDA
72  break;
73 #endif
75  throw memory_exception("not initialised!");
76  default:
77  throw memory_exception("not implemented");
78  }
79 }
80 
86 template<typename NumericT>
89 {
90  switch (viennacl::traits::handle(A).get_active_handle_id())
91  {
94  break;
95 #ifdef VIENNACL_WITH_OPENCL
98  break;
99 #endif
100 #ifdef VIENNACL_WITH_CUDA
103  break;
104 #endif
106  throw memory_exception("not initialised!");
107  default:
108  throw memory_exception("not implemented");
109  }
110 }
111 
120 template<typename NumericT>
122  vector<NumericT> & aij_L)
123 {
124  switch (viennacl::traits::handle(L).get_active_handle_id())
125  {
128  break;
129 #ifdef VIENNACL_WITH_OPENCL
132  break;
133 #endif
134 #ifdef VIENNACL_WITH_CUDA
137  break;
138 #endif
140  throw memory_exception("not initialised!");
141  default:
142  throw memory_exception("not implemented");
143  }
144 }
145 
146 
147 
149 
155 template<typename NumericT>
159 {
160  switch (viennacl::traits::handle(A).get_active_handle_id())
161  {
164  break;
165 #ifdef VIENNACL_WITH_OPENCL
168  break;
169 #endif
170 #ifdef VIENNACL_WITH_CUDA
173  break;
174 #endif
176  throw memory_exception("not initialised!");
177  default:
178  throw memory_exception("not implemented");
179  }
180 }
181 
187 template<typename NumericT>
191 {
192  switch (viennacl::traits::handle(A).get_active_handle_id())
193  {
196  break;
197 #ifdef VIENNACL_WITH_OPENCL
200  break;
201 #endif
202 #ifdef VIENNACL_WITH_CUDA
205  break;
206 #endif
208  throw memory_exception("not initialised!");
209  default:
210  throw memory_exception("not implemented");
211  }
212 }
213 
219 template<typename NumericT>
222 {
225  (void)orig_ctx;
226  (void)cpu_ctx;
227 
228  viennacl::compressed_matrix<NumericT> A_host(0, 0, 0, cpu_ctx);
229  (void)A_host;
230 
232  {
235  break;
236 #ifdef VIENNACL_WITH_OPENCL
238  A_host = A;
239  B.switch_memory_context(cpu_ctx);
241  B.switch_memory_context(orig_ctx);
242  break;
243 #endif
244 #ifdef VIENNACL_WITH_CUDA
246  A_host = A;
247  B.switch_memory_context(cpu_ctx);
249  B.switch_memory_context(orig_ctx);
250  break;
251 #endif
253  throw memory_exception("not initialised!");
254  default:
255  throw memory_exception("not implemented");
256  }
257 }
258 
259 
260 
271 template<typename NumericT>
273  vector<NumericT> const & aij_L,
274  compressed_matrix<NumericT> & U_trans,
275  vector<NumericT> const & aij_U_trans)
276 {
277  switch (viennacl::traits::handle(L).get_active_handle_id())
278  {
280  viennacl::linalg::host_based::ilu_chow_patel_sweep(L, aij_L, U_trans, aij_U_trans);
281  break;
282 #ifdef VIENNACL_WITH_OPENCL
284  viennacl::linalg::opencl::ilu_chow_patel_sweep(L, aij_L, U_trans, aij_U_trans);
285  break;
286 #endif
287 #ifdef VIENNACL_WITH_CUDA
289  viennacl::linalg::cuda::ilu_chow_patel_sweep(L, aij_L, U_trans, aij_U_trans);
290  break;
291 #endif
293  throw memory_exception("not initialised!");
294  default:
295  throw memory_exception("not implemented");
296  }
297 }
298 
304 template<typename NumericT>
306  vector<NumericT> & diag_R)
307 {
308  switch (viennacl::traits::handle(R).get_active_handle_id())
309  {
312  break;
313 #ifdef VIENNACL_WITH_OPENCL
316  break;
317 #endif
318 #ifdef VIENNACL_WITH_CUDA
321  break;
322 #endif
324  throw memory_exception("not initialised!");
325  default:
326  throw memory_exception("not implemented");
327  }
328 }
329 
330 } //namespace linalg
331 } //namespace viennacl
332 
333 
334 #endif
Exception class in case of memory errors.
Definition: forwards.h:572
Generic size and resize functionality for different vector and matrix types.
void extract_LU(compressed_matrix< NumericT > const &A, compressed_matrix< NumericT > &L, compressed_matrix< NumericT > &U)
Extracts the underlying OpenCL start index handle from a vector, a matrix, an expression etc...
Various little tools used here and there in ViennaCL.
void switch_memory_context(viennacl::context new_ctx)
Switches the memory context of the matrix.
void extract_L(compressed_matrix< NumericT > const &A, compressed_matrix< NumericT > &L)
Extracts the lower triangular part L from A.
void ilu_scale(compressed_matrix< NumericT > const &A, compressed_matrix< NumericT > &L, compressed_matrix< NumericT > &U)
Scales the values extracted from A such that A' = DAD has unit diagonal. Updates values from A in L a...
void icc_chow_patel_sweep(compressed_matrix< NumericT > &L, vector< NumericT > &aij_L)
Performs one nonlinear relaxation step in the Chow-Patel-ICC using OpenMP (cf. Algorithm 3 in paper...
void ilu_chow_patel_sweep(compressed_matrix< NumericT > &L, vector< NumericT > const &aij_L, compressed_matrix< NumericT > &U_trans, vector< NumericT > const &aij_U_trans)
Performs one nonlinear relaxation step in the Chow-Patel-ILU (cf. Algorithm 2 in paper) ...
Implementations of specialized routines for the Chow-Patel parallel ILU preconditioner using OpenCL...
This file provides the forward declarations for the main types used within ViennaCL.
Determines row and column increments for matrices and matrix proxies.
void extract_LU(compressed_matrix< NumericT > const &A, compressed_matrix< NumericT > &L, compressed_matrix< NumericT > &U)
void ilu_chow_patel_sweep(compressed_matrix< NumericT > &L, vector< NumericT > const &aij_L, compressed_matrix< NumericT > &U_trans, vector< NumericT > const &aij_U_trans)
Performs one nonlinear relaxation step in the Chow-Patel-ILU using OpenMP (cf. Algorithm 2 in paper) ...
void extract_L(compressed_matrix< NumericT > const &A, compressed_matrix< NumericT > &L)
void ilu_form_neumann_matrix(compressed_matrix< NumericT > &R, vector< NumericT > &diag_R)
Represents a generic 'context' similar to an OpenCL context, but is backend-agnostic and thus also su...
Definition: context.hpp:39
void ilu_scale(compressed_matrix< NumericT > const &A, compressed_matrix< NumericT > &L, compressed_matrix< NumericT > &U)
Scales the values extracted from A such that A' = DAD has unit diagonal. Updates values from A in L a...
Main namespace in ViennaCL. Holds all the basic types such as vector, matrix, etc. and defines operations upon them.
Definition: cpu_ram.hpp:34
void ilu_form_neumann_matrix(compressed_matrix< NumericT > &R, vector< NumericT > &diag_R)
void extract_L(compressed_matrix< NumericT > const &A, compressed_matrix< NumericT > &L)
void extract_L(compressed_matrix< NumericT > const &A, compressed_matrix< NumericT > &L)
void icc_scale(compressed_matrix< NumericT > const &A, compressed_matrix< NumericT > &L)
Scales the values extracted from A such that A' = DAD has unit diagonal. Updates values from A in L a...
void ilu_scale(compressed_matrix< NumericT > const &A, compressed_matrix< NumericT > &L, compressed_matrix< NumericT > &U)
Scales the values extracted from A such that A' = DAD has unit diagonal. Updates values from A in L a...
void extract_LU(compressed_matrix< NumericT > const &A, compressed_matrix< NumericT > &L, compressed_matrix< NumericT > &U)
Extracts the lower triangular part L and the upper triangular part U from A.
void icc_scale(compressed_matrix< NumericT > const &A, compressed_matrix< NumericT > &L)
Scales the values extracted from A such that A' = DAD has unit diagonal. Updates values from A in L a...
void ilu_scale(compressed_matrix< NumericT > const &A, compressed_matrix< NumericT > &L, compressed_matrix< NumericT > &U)
Scales the values extracted from A such that A' = DAD has unit diagonal. Updates values from A in L a...
All the predicates used within ViennaCL. Checks for expressions to be vectors, etc.
void ilu_form_neumann_matrix(compressed_matrix< NumericT > &R, vector< NumericT > &diag_R)
Extracts the lower triangular part L and the upper triangular part U from A.
void ilu_transpose(compressed_matrix< NumericT > const &A, compressed_matrix< NumericT > &B)
Transposition B <- A^T, where the aij-vector is permuted in the same way as the value array in A when...
void icc_scale(compressed_matrix< NumericT > const &A, compressed_matrix< NumericT > &L)
Scales the values extracted from A such that A' = DAD has unit diagonal. Updates values from A in L a...
viennacl::context context(T const &t)
Returns an ID for the currently active memory domain of an object.
Definition: context.hpp:40
void ilu_transpose(compressed_matrix< NumericT > const &A, compressed_matrix< NumericT > &B)
Implementations of specialized routines for the Chow-Patel parallel ILU preconditioner using the host...
void icc_chow_patel_sweep(compressed_matrix< NumericT > &L, vector< NumericT > &aij_L)
Performs one nonlinear relaxation step in the Chow-Patel-ICC (cf. Algorithm 3 in paper, but for L rather than U)
void ilu_form_neumann_matrix(compressed_matrix< NumericT > &R, vector< NumericT > &diag_R)
Implementation of a range object for use with proxy objects.
void icc_chow_patel_sweep(compressed_matrix< NumericT > &L, vector< NumericT > const &aij_L)
Performs one nonlinear relaxation step in the Chow-Patel-ILU using OpenMP (cf. Algorithm 2 in paper) ...
Extracts the underlying OpenCL handle from a vector, a matrix, an expression etc. ...
void ilu_chow_patel_sweep(compressed_matrix< NumericT > &L, vector< NumericT > const &aij_L, compressed_matrix< NumericT > &U_trans, vector< NumericT > const &aij_U_trans)
Performs one nonlinear relaxation step in the Chow-Patel-ILU using OpenMP (cf. Algorithm 2 in paper) ...
viennacl::backend::mem_handle & handle(T &obj)
Returns the generic memory handle of an object. Non-const version.
Definition: handle.hpp:41
void extract_LU(compressed_matrix< NumericT > const &A, compressed_matrix< NumericT > &L, compressed_matrix< NumericT > &U)
Implementations of specialized routines for the Chow-Patel parallel ILU preconditioner using CUDA...
Implementation of the ViennaCL scalar class.
void icc_scale(compressed_matrix< NumericT > const &A, compressed_matrix< NumericT > &L)
Scales the values extracted from A such that A' = DAD has unit diagonal. Updates values from A in L a...
void ilu_chow_patel_sweep(compressed_matrix< NumericT > &L, vector< NumericT > const &aij_L, compressed_matrix< NumericT > &U_trans, vector< NumericT > const &aij_U_trans)
Performs one nonlinear relaxation step in the Chow-Patel-ILU using OpenMP (cf. Algorithm 2 in paper) ...
Simple enable-if variant that uses the SFINAE pattern.
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
void icc_chow_patel_sweep(compressed_matrix< NumericT > &L, vector< NumericT > const &aij_L)
Performs one nonlinear relaxation step in the Chow-Patel-ILU using OpenMP (cf. Algorithm 2 in paper) ...