1 #ifndef VIENNACL_LINALG_BICGSTAB_HPP_
2 #define VIENNACL_LINALG_BICGSTAB_HPP_
57 : tol_(tol), iterations_(max_iters), iterations_before_restart_(max_iters_before_restart) {}
71 double error()
const {
return last_error_; }
73 void error(
double e)
const { last_error_ = e; }
82 mutable double last_error_;
90 template<
typename MatrixT,
typename NumericT>
115 std::vector<NumericT> host_inner_prod_buffer(inner_prod_buffer.
size());
121 NumericT residual_norm = norm_rhs_host;
122 inner_prod_buffer[0] = norm_rhs_host * norm_rhs_host;
124 NumericT r_dot_r0 = 0;
125 NumericT As_dot_As = 0;
126 NumericT As_dot_s = 0;
127 NumericT Ap_dot_r0 = 0;
128 NumericT As_dot_r0 = 0;
129 NumericT s_dot_s = 0;
131 if (norm_rhs_host <= 0)
140 inner_prod_buffer, buffer_size_per_vector, 3*buffer_size_per_vector);
160 inner_prod_buffer, buffer_size_per_vector, 5*buffer_size_per_vector);
167 inner_prod_buffer, buffer_size_per_vector, 4*buffer_size_per_vector);
173 r_dot_r0 = std::accumulate(host_inner_prod_buffer.begin(), host_inner_prod_buffer.begin() + buffer_size_per_vector, NumericT(0));
174 As_dot_As = std::accumulate(host_inner_prod_buffer.begin() + buffer_size_per_vector, host_inner_prod_buffer.begin() + 2 * buffer_size_per_vector, NumericT(0));
175 As_dot_s = std::accumulate(host_inner_prod_buffer.begin() + 2 * buffer_size_per_vector, host_inner_prod_buffer.begin() + 3 * buffer_size_per_vector, NumericT(0));
176 Ap_dot_r0 = std::accumulate(host_inner_prod_buffer.begin() + 3 * buffer_size_per_vector, host_inner_prod_buffer.begin() + 4 * buffer_size_per_vector, NumericT(0));
177 As_dot_r0 = std::accumulate(host_inner_prod_buffer.begin() + 4 * buffer_size_per_vector, host_inner_prod_buffer.begin() + 5 * buffer_size_per_vector, NumericT(0));
178 s_dot_s = std::accumulate(host_inner_prod_buffer.begin() + 5 * buffer_size_per_vector, host_inner_prod_buffer.begin() + 6 * buffer_size_per_vector, NumericT(0));
180 alpha = r_dot_r0 / Ap_dot_r0;
181 beta = - As_dot_r0 / Ap_dot_r0;
182 omega = As_dot_s / As_dot_As;
184 residual_norm = std::sqrt(s_dot_s - NumericT(2.0) * omega * As_dot_s + omega * omega * As_dot_As);
185 if (std::fabs(residual_norm / norm_rhs_host) < tag.
tolerance())
195 r0star, inner_prod_buffer, buffer_size_per_vector);
199 tag.
error(residual_norm / norm_rhs_host);
208 template<
typename NumericT>
218 template<
typename NumericT>
228 template<
typename NumericT>
238 template<
typename NumericT>
251 template<
typename NumericT>
261 template<
typename NumericT>
271 template<
typename NumericT>
281 template<
typename NumericT>
293 template<
typename NumericT>
303 template<
typename NumericT>
313 template<
typename NumericT>
323 template<
typename NumericT>
335 template<
typename NumericT>
345 template<
typename NumericT>
355 template<
typename NumericT>
365 template<
typename NumericT>
378 template<
typename NumericT>
388 template<
typename NumericT>
398 template<
typename NumericT>
408 template<
typename NumericT>
430 template<
typename MatrixT,
typename VectorT>
435 VectorT result = rhs;
438 VectorT residual = rhs;
440 VectorT r0star = rhs;
446 CPU_NumericType ip_rr0star = norm_rhs_host * norm_rhs_host;
447 CPU_NumericType beta;
448 CPU_NumericType alpha;
449 CPU_NumericType omega;
451 CPU_NumericType new_ip_rr0star = 0;
452 CPU_NumericType residual_norm = norm_rhs_host;
454 if (norm_rhs_host <= 0)
457 bool restart_flag =
true;
468 ip_rr0star *= ip_rr0star;
469 restart_flag =
false;
477 s = residual - alpha*tmp0;
483 result += alpha * p + omega * s;
484 residual = s - omega * tmp1;
488 if (std::fabs(residual_norm / norm_rhs_host) < tag.
tolerance())
491 beta = new_ip_rr0star / ip_rr0star * alpha/omega;
492 ip_rr0star = new_ip_rr0star;
494 if ( (ip_rr0star <= 0 && ip_rr0star >= 0)
495 || (omega <= 0 && omega >= 0)
504 p = residual + beta * p;
508 tag.
error(residual_norm / norm_rhs_host);
513 template<
typename MatrixT,
typename VectorT>
516 return solve(matrix, rhs, tag);
529 template<
typename MatrixT,
typename VectorT,
typename PreconditionerT>
534 VectorT result = rhs;
537 VectorT residual = rhs;
538 VectorT r0star = residual;
543 VectorT p = residual;
547 CPU_NumericType beta;
548 CPU_NumericType alpha;
549 CPU_NumericType omega;
550 CPU_NumericType new_ip_rr0star = 0;
551 CPU_NumericType residual_norm = norm_rhs_host;
556 bool restart_flag =
true;
564 precond.apply(residual);
568 ip_rr0star *= ip_rr0star;
569 restart_flag =
false;
578 s = residual - alpha*tmp0;
585 result += alpha * p + omega * s;
586 residual = s - omega * tmp1;
589 if (residual_norm / norm_rhs_host < tag.
tolerance())
594 beta = new_ip_rr0star / ip_rr0star * alpha/omega;
595 ip_rr0star = new_ip_rr0star;
604 p = residual + beta * p;
610 tag.
error(residual_norm / norm_rhs_host);
Sparse matrix class using a hybrid format composed of the ELL and CSR format for storing the nonzeros...
vcl_size_t max_iterations_before_restart() const
Returns the maximum number of iterations before a restart.
T norm_2(std::vector< T, A > const &v1)
viennacl::vector< NumericT > pipelined_solve(MatrixT const &A, viennacl::vector_base< NumericT > const &rhs, bicgstab_tag const &tag, viennacl::linalg::no_precond)
Implementation of a pipelined stabilized Bi-conjugate gradient solver.
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.
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...
double error() const
Returns the estimated relative error at the end of the solver run.
void clear(VectorType &vec)
Generic routine for setting all entries of a vector to zero. This is the version for non-ViennaCL obj...
This file provides the forward declarations for the main types used within ViennaCL.
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)
Generic interface for the computation of inner products. See viennacl/linalg/vector_operations.hpp for implementations.
double tolerance() const
Returns the relative tolerance.
viennacl::vector< NumericT > solve(viennacl::compressed_matrix< NumericT > const &A, viennacl::vector_base< NumericT > const &rhs, bicgstab_tag const &tag, viennacl::linalg::no_precond)
Overload for the pipelined BiCGStab implementation for the ViennaCL sparse matrix types...
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)
Class for representing non-strided subvectors of a bigger vector x.
Sparse matrix class using the ELLPACK format for storing the nonzeros.
iterator begin()
Returns an iterator pointing to the beginning of the vector (STL like)
A tag class representing the use of no preconditioner.
Sparse matrix class using the sliced ELLPACK with parameters C, .
Extracts the underlying context from objects.
Class for representing strided subvectors of a bigger vector x.
vcl_size_t max_iterations() const
Returns the maximum number of iterations.
Generic clear functionality for different vector and matrix types.
T::ERROR_CANNOT_DEDUCE_CPU_SCALAR_TYPE_FOR_T type
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.
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.
void error(double e) const
Sets the estimated relative error at the end of the solver run.
size_type size() const
Returns the length of the vector (cf. std::vector)
vcl_size_t iters() const
Return the number of solver iterations:
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
A tag for the stabilized Bi-conjugate gradient solver. Used for supplying solver parameters and for d...
iterator end()
Returns an iterator pointing to the end of the vector (STL like)
A collection of compile time type deductions.
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)