ViennaCL - The Vienna Computing Library  1.7.0
Free open-source GPU-accelerated linear algebra and solver library.
bicgstab.hpp
Go to the documentation of this file.
1 #ifndef VIENNACL_LINALG_BICGSTAB_HPP_
2 #define VIENNACL_LINALG_BICGSTAB_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 <vector>
26 #include <cmath>
27 #include <numeric>
28 
29 #include "viennacl/forwards.h"
30 #include "viennacl/tools/tools.hpp"
31 #include "viennacl/linalg/prod.hpp"
35 #include "viennacl/traits/size.hpp"
39 
40 namespace viennacl
41 {
42 namespace linalg
43 {
44 
48 {
49 public:
56  bicgstab_tag(double tol = 1e-8, vcl_size_t max_iters = 400, vcl_size_t max_iters_before_restart = 200)
57  : tol_(tol), abs_tol_(0), iterations_(max_iters), iterations_before_restart_(max_iters_before_restart) {}
58 
60  double tolerance() const { return tol_; }
61 
63  double abs_tolerance() const { return abs_tol_; }
65  void abs_tolerance(double new_tol) { if (new_tol >= 0) abs_tol_ = new_tol; }
66 
68  vcl_size_t max_iterations() const { return iterations_; }
70  vcl_size_t max_iterations_before_restart() const { return iterations_before_restart_; }
71 
73  vcl_size_t iters() const { return iters_taken_; }
74  void iters(vcl_size_t i) const { iters_taken_ = i; }
75 
77  double error() const { return last_error_; }
79  void error(double e) const { last_error_ = e; }
80 
81 private:
82  double tol_;
83  double abs_tol_;
84  vcl_size_t iterations_;
85  vcl_size_t iterations_before_restart_;
86 
87  //return values from solver
88  mutable vcl_size_t iters_taken_;
89  mutable double last_error_;
90 };
91 
92 
93 
94 namespace detail
95 {
97  template<typename MatrixT, typename NumericT>
98  viennacl::vector<NumericT> pipelined_solve(MatrixT const & A, //MatrixType const & A,
100  bicgstab_tag const & tag,
102  bool (*monitor)(viennacl::vector<NumericT> const &, NumericT, void*) = NULL,
103  void *monitor_data = NULL)
104  {
106 
107  viennacl::vector<NumericT> residual = rhs;
109  viennacl::vector<NumericT> r0star = rhs;
113 
114  // Layout of temporary buffer:
115  // chunk 0: <residual, r_0^*>
116  // chunk 1: <As, As>
117  // chunk 2: <As, s>
118  // chunk 3: <Ap, r_0^*>
119  // chunk 4: <As, r_0^*>
120  // chunk 5: <s, s>
121  vcl_size_t buffer_size_per_vector = 256;
122  vcl_size_t num_buffer_chunks = 6;
123  viennacl::vector<NumericT> inner_prod_buffer = viennacl::zero_vector<NumericT>(num_buffer_chunks*buffer_size_per_vector, viennacl::traits::context(rhs)); // temporary buffer
124  std::vector<NumericT> host_inner_prod_buffer(inner_prod_buffer.size());
125 
126  NumericT norm_rhs_host = viennacl::linalg::norm_2(residual);
127  NumericT beta;
128  NumericT alpha;
129  NumericT omega;
130  NumericT residual_norm = norm_rhs_host;
131  inner_prod_buffer[0] = norm_rhs_host * norm_rhs_host;
132 
133  NumericT r_dot_r0 = 0;
134  NumericT As_dot_As = 0;
135  NumericT As_dot_s = 0;
136  NumericT Ap_dot_r0 = 0;
137  NumericT As_dot_r0 = 0;
138  NumericT s_dot_s = 0;
139 
140  if (norm_rhs_host <= tag.abs_tolerance()) //solution is zero if RHS norm is zero
141  return result;
142 
143  for (vcl_size_t i = 0; i < tag.max_iterations(); ++i)
144  {
145  tag.iters(i+1);
146  // Ap = A*p_j
147  // Ap_dot_r0 = <Ap, r_0^*>
149  inner_prod_buffer, buffer_size_per_vector, 3*buffer_size_per_vector);
150 
152 
154  //
156  //viennacl::fast_copy(inner_prod_buffer.begin(), inner_prod_buffer.end(), host_inner_prod_buffer.begin());
157  //Ap_dot_r0 = std::accumulate(host_inner_prod_buffer.begin() + buffer_size_per_vector, host_inner_prod_buffer.begin() + 2 * buffer_size_per_vector, ScalarType(0));
158 
159  //alpha = residual_dot_r0 / Ap_dot_r0;
160 
162  //s = residual - alpha * Ap;
163 
165  // s = r - alpha * Ap
166  // <s, s> first stage
167  // dump alpha at end of inner_prod_buffer
169  inner_prod_buffer, buffer_size_per_vector, 5*buffer_size_per_vector);
170 
171  // As = A*s_j
172  // As_dot_As = <As, As>
173  // As_dot_s = <As, s>
174  // As_dot_r0 = <As, r_0^*>
176  inner_prod_buffer, buffer_size_per_vector, 4*buffer_size_per_vector);
177 
179 
180  viennacl::fast_copy(inner_prod_buffer.begin(), inner_prod_buffer.end(), host_inner_prod_buffer.begin());
181 
182  typedef typename std::vector<NumericT>::difference_type difference_type;
183 
184  r_dot_r0 = std::accumulate(host_inner_prod_buffer.begin(), host_inner_prod_buffer.begin() + difference_type( buffer_size_per_vector), NumericT(0));
185  As_dot_As = std::accumulate(host_inner_prod_buffer.begin() + difference_type( buffer_size_per_vector), host_inner_prod_buffer.begin() + difference_type(2 * buffer_size_per_vector), NumericT(0));
186  As_dot_s = std::accumulate(host_inner_prod_buffer.begin() + difference_type(2 * buffer_size_per_vector), host_inner_prod_buffer.begin() + difference_type(3 * buffer_size_per_vector), NumericT(0));
187  Ap_dot_r0 = std::accumulate(host_inner_prod_buffer.begin() + difference_type(3 * buffer_size_per_vector), host_inner_prod_buffer.begin() + difference_type(4 * buffer_size_per_vector), NumericT(0));
188  As_dot_r0 = std::accumulate(host_inner_prod_buffer.begin() + difference_type(4 * buffer_size_per_vector), host_inner_prod_buffer.begin() + difference_type(5 * buffer_size_per_vector), NumericT(0));
189  s_dot_s = std::accumulate(host_inner_prod_buffer.begin() + difference_type(5 * buffer_size_per_vector), host_inner_prod_buffer.begin() + difference_type(6 * buffer_size_per_vector), NumericT(0));
190 
191  alpha = r_dot_r0 / Ap_dot_r0;
192  beta = - As_dot_r0 / Ap_dot_r0;
193  omega = As_dot_s / As_dot_As;
194 
195  residual_norm = std::sqrt(s_dot_s - NumericT(2.0) * omega * As_dot_s + omega * omega * As_dot_As);
196  if (monitor && monitor(result, std::fabs(residual_norm / norm_rhs_host), monitor_data))
197  break;
198  if (std::fabs(residual_norm / norm_rhs_host) < tag.tolerance() || residual_norm < tag.abs_tolerance())
199  break;
200 
201  // x_{j+1} = x_j + alpha * p_j + omega * s_j
202  // r_{j+1} = s_j - omega * t_j
203  // p_{j+1} = r_{j+1} + beta * (p_j - omega * q_j)
204  // and compute first stage of r_dot_r0 = <r_{j+1}, r_o^*> for use in next iteration
205  viennacl::linalg::pipelined_bicgstab_vector_update(result, alpha, p, omega, s,
206  residual, As,
207  beta, Ap,
208  r0star, inner_prod_buffer, buffer_size_per_vector);
209  }
210 
211  //store last error estimate:
212  tag.error(residual_norm / norm_rhs_host);
213 
214  return result;
215  }
216 
218  template<typename NumericT>
220  viennacl::vector<NumericT> const & rhs,
221  bicgstab_tag const & tag,
223  bool (*monitor)(viennacl::vector<NumericT> const &, NumericT, void*) = NULL,
224  void *monitor_data = NULL)
225  {
226  return pipelined_solve(A, rhs, tag, viennacl::linalg::no_precond(), monitor, monitor_data);
227  }
228 
229 
231  template<typename NumericT>
233  viennacl::vector<NumericT> const & rhs,
234  bicgstab_tag const & tag,
236  bool (*monitor)(viennacl::vector<NumericT> const &, NumericT, void*) = NULL,
237  void *monitor_data = NULL)
238  {
239  return detail::pipelined_solve(A, rhs, tag, viennacl::linalg::no_precond(), monitor, monitor_data);
240  }
241 
242 
243 
245  template<typename NumericT>
247  viennacl::vector<NumericT> const & rhs,
248  bicgstab_tag const & tag,
250  bool (*monitor)(viennacl::vector<NumericT> const &, NumericT, void*) = NULL,
251  void *monitor_data = NULL)
252  {
253  return detail::pipelined_solve(A, rhs, tag, viennacl::linalg::no_precond(), monitor, monitor_data);
254  }
255 
256 
257 
259  template<typename NumericT>
261  viennacl::vector<NumericT> const & rhs,
262  bicgstab_tag const & tag,
264  bool (*monitor)(viennacl::vector<NumericT> const &, NumericT, void*) = NULL,
265  void *monitor_data = NULL)
266  {
267  return detail::pipelined_solve(A, rhs, tag, viennacl::linalg::no_precond(), monitor, monitor_data);
268  }
269 
270 
272  template<typename NumericT>
274  viennacl::vector<NumericT> const & rhs,
275  bicgstab_tag const & tag,
277  bool (*monitor)(viennacl::vector<NumericT> const &, NumericT, void*) = NULL,
278  void *monitor_data = NULL)
279  {
280  return detail::pipelined_solve(A, rhs, tag, viennacl::linalg::no_precond(), monitor, monitor_data);
281  }
282 
283 
295  template<typename MatrixT, typename VectorT>
296  VectorT solve_impl(MatrixT const & matrix,
297  VectorT const & rhs,
298  bicgstab_tag const & tag,
300  bool (*monitor)(VectorT const &, typename viennacl::result_of::cpu_value_type<typename viennacl::result_of::value_type<VectorT>::type>::type, void*) = NULL,
301  void *monitor_data = NULL)
302  {
303  typedef typename viennacl::result_of::value_type<VectorT>::type NumericType;
304  typedef typename viennacl::result_of::cpu_value_type<NumericType>::type CPU_NumericType;
305  VectorT result = rhs;
306  viennacl::traits::clear(result);
307 
308  VectorT residual = rhs;
309  VectorT p = rhs;
310  VectorT r0star = rhs;
311  VectorT tmp0 = rhs;
312  VectorT tmp1 = rhs;
313  VectorT s = rhs;
314 
315  CPU_NumericType norm_rhs_host = viennacl::linalg::norm_2(residual);
316  CPU_NumericType ip_rr0star = norm_rhs_host * norm_rhs_host;
317  CPU_NumericType beta;
318  CPU_NumericType alpha;
319  CPU_NumericType omega;
320  //ScalarType inner_prod_temp; //temporary variable for inner product computation
321  CPU_NumericType new_ip_rr0star = 0;
322  CPU_NumericType residual_norm = norm_rhs_host;
323 
324  if (norm_rhs_host <= tag.abs_tolerance()) //solution is zero if RHS norm is zero
325  return result;
326 
327  bool restart_flag = true;
328  vcl_size_t last_restart = 0;
329  for (vcl_size_t i = 0; i < tag.max_iterations(); ++i)
330  {
331  if (restart_flag)
332  {
333  residual = viennacl::linalg::prod(matrix, result);
334  residual = rhs - residual;
335  p = residual;
336  r0star = residual;
337  ip_rr0star = viennacl::linalg::norm_2(residual);
338  ip_rr0star *= ip_rr0star;
339  restart_flag = false;
340  last_restart = i;
341  }
342 
343  tag.iters(i+1);
344  tmp0 = viennacl::linalg::prod(matrix, p);
345  alpha = ip_rr0star / viennacl::linalg::inner_prod(tmp0, r0star);
346 
347  s = residual - alpha*tmp0;
348 
349  tmp1 = viennacl::linalg::prod(matrix, s);
350  CPU_NumericType norm_tmp1 = viennacl::linalg::norm_2(tmp1);
351  omega = viennacl::linalg::inner_prod(tmp1, s) / (norm_tmp1 * norm_tmp1);
352 
353  result += alpha * p + omega * s;
354  residual = s - omega * tmp1;
355 
356  new_ip_rr0star = viennacl::linalg::inner_prod(residual, r0star);
357  residual_norm = viennacl::linalg::norm_2(residual);
358  if (monitor && monitor(result, std::fabs(residual_norm / norm_rhs_host), monitor_data))
359  break;
360  if (std::fabs(residual_norm / norm_rhs_host) < tag.tolerance() || residual_norm < tag.abs_tolerance())
361  break;
362 
363  beta = new_ip_rr0star / ip_rr0star * alpha/omega;
364  ip_rr0star = new_ip_rr0star;
365 
366  if ( (ip_rr0star <= 0 && ip_rr0star >= 0)
367  || (omega <= 0 && omega >= 0)
368  || (i - last_restart > tag.max_iterations_before_restart())
369  ) //search direction degenerate. A restart might help
370  restart_flag = true;
371 
372  // Execution of
373  // p = residual + beta * (p - omega*tmp0);
374  // without introducing temporary vectors:
375  p -= omega * tmp0;
376  p = residual + beta * p;
377  }
378 
379  //store last error estimate:
380  tag.error(residual_norm / norm_rhs_host);
381 
382  return result;
383  }
384 
385 
398  template<typename MatrixT, typename VectorT, typename PreconditionerT>
399  VectorT solve_impl(MatrixT const & matrix,
400  VectorT const & rhs,
401  bicgstab_tag const & tag,
402  PreconditionerT const & precond,
403  bool (*monitor)(VectorT const &, typename viennacl::result_of::cpu_value_type<typename viennacl::result_of::value_type<VectorT>::type>::type, void*) = NULL,
404  void *monitor_data = NULL)
405  {
406  typedef typename viennacl::result_of::value_type<VectorT>::type NumericType;
407  typedef typename viennacl::result_of::cpu_value_type<NumericType>::type CPU_NumericType;
408  VectorT result = rhs;
409  viennacl::traits::clear(result);
410 
411  VectorT residual = rhs;
412  VectorT r0star = residual; //can be chosen arbitrarily in fact
413  VectorT tmp0 = rhs;
414  VectorT tmp1 = rhs;
415  VectorT s = rhs;
416 
417  VectorT p = residual;
418 
419  CPU_NumericType ip_rr0star = viennacl::linalg::norm_2(residual);
420  CPU_NumericType norm_rhs_host = viennacl::linalg::norm_2(residual);
421  CPU_NumericType beta;
422  CPU_NumericType alpha;
423  CPU_NumericType omega;
424  CPU_NumericType new_ip_rr0star = 0;
425  CPU_NumericType residual_norm = norm_rhs_host;
426 
427  if (norm_rhs_host <= tag.abs_tolerance()) //solution is zero if RHS norm is zero
428  return result;
429 
430  bool restart_flag = true;
431  vcl_size_t last_restart = 0;
432  for (unsigned int i = 0; i < tag.max_iterations(); ++i)
433  {
434  if (restart_flag)
435  {
436  residual = viennacl::linalg::prod(matrix, result);
437  residual = rhs - residual;
438  precond.apply(residual);
439  p = residual;
440  r0star = residual;
441  ip_rr0star = viennacl::linalg::norm_2(residual);
442  ip_rr0star *= ip_rr0star;
443  restart_flag = false;
444  last_restart = i;
445  }
446 
447  tag.iters(i+1);
448  tmp0 = viennacl::linalg::prod(matrix, p);
449  precond.apply(tmp0);
450  alpha = ip_rr0star / viennacl::linalg::inner_prod(tmp0, r0star);
451 
452  s = residual - alpha*tmp0;
453 
454  tmp1 = viennacl::linalg::prod(matrix, s);
455  precond.apply(tmp1);
456  CPU_NumericType norm_tmp1 = viennacl::linalg::norm_2(tmp1);
457  omega = viennacl::linalg::inner_prod(tmp1, s) / (norm_tmp1 * norm_tmp1);
458 
459  result += alpha * p + omega * s;
460  residual = s - omega * tmp1;
461 
462  residual_norm = viennacl::linalg::norm_2(residual);
463  if (monitor && monitor(result, std::fabs(residual_norm / norm_rhs_host), monitor_data))
464  break;
465  if (residual_norm / norm_rhs_host < tag.tolerance() || residual_norm < tag.abs_tolerance())
466  break;
467 
468  new_ip_rr0star = viennacl::linalg::inner_prod(residual, r0star);
469 
470  beta = new_ip_rr0star / ip_rr0star * alpha/omega;
471  ip_rr0star = new_ip_rr0star;
472 
473  if (!ip_rr0star || !omega || i - last_restart > tag.max_iterations_before_restart()) //search direction degenerate. A restart might help
474  restart_flag = true;
475 
476  // Execution of
477  // p = residual + beta * (p - omega*tmp0);
478  // without introducing temporary vectors:
479  p -= omega * tmp0;
480  p = residual + beta * p;
481 
482  //std::cout << "Rel. Residual in current step: " << std::sqrt(std::fabs(viennacl::linalg::inner_prod(residual, residual) / norm_rhs_host)) << std::endl;
483  }
484 
485  //store last error estimate:
486  tag.error(residual_norm / norm_rhs_host);
487 
488  return result;
489  }
490 
491 }
492 
493 
494 
495 template<typename MatrixT, typename VectorT, typename PreconditionerT>
496 VectorT solve(MatrixT const & matrix, VectorT const & rhs, bicgstab_tag const & tag, PreconditionerT const & precond)
497 {
498  return detail::solve_impl(matrix, rhs, tag, precond);
499 }
500 
501 
507 template<typename IndexT, typename NumericT, typename PreconditionerT>
508 std::vector<NumericT> solve(std::vector< std::map<IndexT, NumericT> > const & A, std::vector<NumericT> const & rhs, bicgstab_tag const & tag, PreconditionerT const & precond)
509 {
511  viennacl::copy(A, vcl_A);
512 
513  viennacl::vector<NumericT> vcl_rhs(rhs.size());
514  viennacl::copy(rhs, vcl_rhs);
515 
516  viennacl::vector<NumericT> vcl_result = solve(vcl_A, vcl_rhs, tag, precond);
517 
518  std::vector<NumericT> result(vcl_result.size());
519  viennacl::copy(vcl_result, result);
520  return result;
521 }
522 
529 template<typename MatrixT, typename VectorT>
530 VectorT solve(MatrixT const & matrix, VectorT const & rhs, bicgstab_tag const & tag)
531 {
532  return solve(matrix, rhs, tag, viennacl::linalg::no_precond());
533 }
534 
535 
536 
537 template<typename VectorT>
539 {
540 public:
542 
543  bicgstab_solver(bicgstab_tag const & tag) : tag_(tag), monitor_callback_(NULL), user_data_(NULL) {}
544 
545  template<typename MatrixT, typename PreconditionerT>
546  VectorT operator()(MatrixT const & A, VectorT const & b, PreconditionerT const & precond) const
547  {
548  if (viennacl::traits::size(init_guess_) > 0) // take initial guess into account
549  {
550  VectorT mod_rhs = viennacl::linalg::prod(A, init_guess_);
551  mod_rhs = b - mod_rhs;
552  VectorT y = detail::solve_impl(A, mod_rhs, tag_, precond, monitor_callback_, user_data_);
553  return init_guess_ + y;
554  }
555  return detail::solve_impl(A, b, tag_, precond, monitor_callback_, user_data_);
556  }
557 
558 
559  template<typename MatrixT>
560  VectorT operator()(MatrixT const & A, VectorT const & b) const
561  {
563  }
564 
569  void set_initial_guess(VectorT const & x) { init_guess_ = x; }
570 
578  void set_monitor(bool (*monitor_fun)(VectorT const &, numeric_type, void *), void *user_data)
579  {
580  monitor_callback_ = monitor_fun;
581  user_data_ = user_data;
582  }
583 
585  bicgstab_tag const & tag() const { return tag_; }
586 
587 private:
588  bicgstab_tag tag_;
589  VectorT init_guess_;
590  bool (*monitor_callback_)(VectorT const &, numeric_type, void *);
591  void *user_data_;
592 };
593 
594 
595 }
596 }
597 
598 #endif
Sparse matrix class using a hybrid format composed of the ELL and CSR format for storing the nonzeros...
Definition: forwards.h:406
vcl_size_t max_iterations_before_restart() const
Returns the maximum number of iterations before a restart.
Definition: bicgstab.hpp:70
T norm_2(std::vector< T, A > const &v1)
Definition: norm_2.hpp:96
Generic interface for the l^2-norm. See viennacl/linalg/vector_operations.hpp for implementations...
bicgstab_tag(double tol=1e-8, vcl_size_t max_iters=400, vcl_size_t max_iters_before_restart=200)
The constructor.
Definition: bicgstab.hpp:56
void abs_tolerance(double new_tol)
Sets the absolute tolerance.
Definition: bicgstab.hpp:65
Generic size and resize functionality for different vector and matrix types.
Generic interface for matrix-vector and matrix-matrix products. See viennacl/linalg/vector_operations...
Various little tools used here and there in ViennaCL.
viennacl::vector< NumericT > solve_impl(viennacl::compressed_matrix< NumericT > const &A, viennacl::vector< NumericT > const &rhs, bicgstab_tag const &tag, viennacl::linalg::no_precond, bool(*monitor)(viennacl::vector< NumericT > const &, NumericT, void *)=NULL, void *monitor_data=NULL)
Overload for the pipelined CG implementation for the ViennaCL sparse matrix types.
Definition: bicgstab.hpp:219
double error() const
Returns the estimated relative error at the end of the solver run.
Definition: bicgstab.hpp:77
void clear(VectorType &vec)
Generic routine for setting all entries of a vector to zero. This is the version for non-ViennaCL obj...
Definition: clear.hpp:43
This file provides the forward declarations for the main types used within ViennaCL.
A dense matrix class.
Definition: forwards.h:375
viennacl::enable_if< viennacl::is_stl< typename viennacl::traits::tag_of< VectorT1 >::type >::value, typename VectorT1::value_type >::type inner_prod(VectorT1 const &v1, VectorT2 const &v2)
Definition: inner_prod.hpp:100
bicgstab_solver(bicgstab_tag const &tag)
Definition: bicgstab.hpp:543
double abs_tolerance() const
Returns the absolute tolerance.
Definition: bicgstab.hpp:63
Generic interface for the computation of inner products. See viennacl/linalg/vector_operations.hpp for implementations.
double tolerance() const
Returns the relative tolerance.
Definition: bicgstab.hpp:60
VectorT solve(MatrixT const &matrix, VectorT const &rhs, bicgstab_tag const &tag, PreconditionerT const &precond)
Definition: bicgstab.hpp:496
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
void pipelined_bicgstab_vector_update(vector_base< NumericT > &result, NumericT alpha, vector_base< NumericT > &p, NumericT omega, vector_base< NumericT > const &s, vector_base< NumericT > &residual, vector_base< NumericT > const &As, NumericT beta, vector_base< NumericT > const &Ap, vector_base< NumericT > const &r0star, vector_base< NumericT > &inner_prod_buffer, vcl_size_t buffer_chunk_size)
Performs a joint vector update operation needed for an efficient pipelined BiCGStab algorithm...
VectorT prod(std::vector< std::vector< T, A1 >, A2 > const &matrix, VectorT const &vector)
Definition: prod.hpp:102
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
Sparse matrix class using the ELLPACK format for storing the nonzeros.
Definition: ell_matrix.hpp:53
iterator begin()
Returns an iterator pointing to the beginning of the vector (STL like)
Definition: vector.hpp:841
A tag class representing the use of no preconditioner.
Definition: forwards.h:873
VectorT operator()(MatrixT const &A, VectorT const &b) const
Definition: bicgstab.hpp:560
Sparse matrix class using the sliced ELLPACK with parameters C, .
Definition: forwards.h:403
Extracts the underlying context from objects.
void set_monitor(bool(*monitor_fun)(VectorT const &, numeric_type, void *), void *user_data)
Sets a monitor function pointer to be called in each iteration. Set to NULL to run without monitor...
Definition: bicgstab.hpp:578
bicgstab_tag const & tag() const
Returns the solver tag containing basic configuration such as tolerances, etc.
Definition: bicgstab.hpp:585
viennacl::vector< NumericT > pipelined_solve(MatrixT const &A, viennacl::vector_base< NumericT > const &rhs, bicgstab_tag const &tag, viennacl::linalg::no_precond, bool(*monitor)(viennacl::vector< NumericT > const &, NumericT, void *)=NULL, void *monitor_data=NULL)
Implementation of a pipelined stabilized Bi-conjugate gradient solver.
Definition: bicgstab.hpp:98
vcl_size_t max_iterations() const
Returns the maximum number of iterations.
Definition: bicgstab.hpp:68
std::size_t vcl_size_t
Definition: forwards.h:75
Generic clear functionality for different vector and matrix types.
T::ERROR_CANNOT_DEDUCE_CPU_SCALAR_TYPE_FOR_T type
Definition: result_of.hpp:271
void pipelined_bicgstab_prod(MatrixT const &A, vector_base< NumericT > const &p, vector_base< NumericT > &Ap, vector_base< NumericT > const &r0star, vector_base< NumericT > &inner_prod_buffer, vcl_size_t buffer_chunk_size, vcl_size_t buffer_chunk_offset)
Performs a joint vector update operation needed for an efficient pipelined CG algorithm.
void set_initial_guess(VectorT const &x)
Specifies an initial guess for the iterative solver.
Definition: bicgstab.hpp:569
Implementations of specialized routines for the iterative solvers.
viennacl::context context(T const &t)
Returns an ID for the currently active memory domain of an object.
Definition: context.hpp:40
viennacl::result_of::cpu_value_type< VectorT >::type numeric_type
Definition: bicgstab.hpp:541
Helper meta function for retrieving the main RAM-based value type. Particularly important to obtain T...
Definition: result_of.hpp:269
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) ...
void error(double e) const
Sets the estimated relative error at the end of the solver run.
Definition: bicgstab.hpp:79
size_type size() const
Returns the length of the vector (cf. std::vector)
Definition: vector_def.hpp:118
vcl_size_t iters() const
Return the number of solver iterations:
Definition: bicgstab.hpp:73
void pipelined_bicgstab_update_s(vector_base< NumericT > &s, vector_base< NumericT > &r, vector_base< NumericT > const &Ap, vector_base< NumericT > &inner_prod_buffer, vcl_size_t buffer_chunk_size, vcl_size_t buffer_chunk_offset)
Performs a joint vector update operation needed for an efficient pipelined CG algorithm.
void iters(vcl_size_t i) const
Definition: bicgstab.hpp:74
A tag for the stabilized Bi-conjugate gradient solver. Used for supplying solver parameters and for d...
Definition: bicgstab.hpp:47
iterator end()
Returns an iterator pointing to the end of the vector (STL like)
Definition: vector.hpp:848
A collection of compile time type deductions.
VectorT operator()(MatrixT const &A, VectorT const &b, PreconditionerT const &precond) const
Definition: bicgstab.hpp:546
A sparse square matrix, where entries are stored as triplets (i,j, val), where i and j are the row an...
void fast_copy(const const_vector_iterator< SCALARTYPE, ALIGNMENT > &gpu_begin, const const_vector_iterator< SCALARTYPE, ALIGNMENT > &gpu_end, CPU_ITERATOR cpu_begin)