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_;
87 template<
typename MatrixT,
typename NumericT>
112 std::vector<NumericT> host_inner_prod_buffer(inner_prod_buffer.
size());
118 NumericT residual_norm = norm_rhs_host;
119 inner_prod_buffer[0] = norm_rhs_host * norm_rhs_host;
121 NumericT r_dot_r0 = 0;
122 NumericT As_dot_As = 0;
123 NumericT As_dot_s = 0;
124 NumericT Ap_dot_r0 = 0;
125 NumericT As_dot_r0 = 0;
126 NumericT s_dot_s = 0;
128 if (norm_rhs_host <= 0)
137 inner_prod_buffer, buffer_size_per_vector, 3*buffer_size_per_vector);
157 inner_prod_buffer, buffer_size_per_vector, 5*buffer_size_per_vector);
164 inner_prod_buffer, buffer_size_per_vector, 4*buffer_size_per_vector);
170 r_dot_r0 = std::accumulate(host_inner_prod_buffer.begin(), host_inner_prod_buffer.begin() + buffer_size_per_vector, NumericT(0));
171 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));
172 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));
173 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));
174 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));
175 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));
177 alpha = r_dot_r0 / Ap_dot_r0;
178 beta = -1.0 * As_dot_r0 / Ap_dot_r0;
179 omega = As_dot_s / As_dot_As;
181 residual_norm = std::sqrt(s_dot_s - 2.0 * omega * As_dot_s + omega * omega * As_dot_As);
182 if (std::fabs(residual_norm / norm_rhs_host) < tag.
tolerance())
192 r0star, inner_prod_buffer, buffer_size_per_vector);
196 tag.
error(residual_norm / norm_rhs_host);
211 template<
typename MatrixT,
typename VectorT>
216 VectorT result = rhs;
219 VectorT residual = rhs;
221 VectorT r0star = rhs;
227 CPU_NumericType ip_rr0star = norm_rhs_host * norm_rhs_host;
228 CPU_NumericType beta;
229 CPU_NumericType alpha;
230 CPU_NumericType omega;
232 CPU_NumericType new_ip_rr0star = 0;
233 CPU_NumericType residual_norm = norm_rhs_host;
235 if (norm_rhs_host <= 0)
238 bool restart_flag =
true;
249 ip_rr0star *= ip_rr0star;
250 restart_flag =
false;
258 s = residual - alpha*tmp0;
264 result += alpha * p + omega * s;
265 residual = s - omega * tmp1;
269 if (std::fabs(residual_norm / norm_rhs_host) < tag.
tolerance())
272 beta = new_ip_rr0star / ip_rr0star * alpha/omega;
273 ip_rr0star = new_ip_rr0star;
275 if ( (ip_rr0star <= 0 && ip_rr0star >= 0)
276 || (omega <= 0 && omega >= 0)
285 p = residual + beta * p;
289 tag.
error(residual_norm / norm_rhs_host);
294 template<
typename MatrixT,
typename VectorT>
297 return solve(matrix, rhs, tag);
310 template<
typename MatrixT,
typename VectorT,
typename PreconditionerT>
315 VectorT result = rhs;
318 VectorT residual = rhs;
319 VectorT r0star = residual;
324 VectorT p = residual;
328 CPU_NumericType beta;
329 CPU_NumericType alpha;
330 CPU_NumericType omega;
331 CPU_NumericType new_ip_rr0star = 0;
332 CPU_NumericType residual_norm = norm_rhs_host;
337 bool restart_flag =
true;
345 precond.apply(residual);
349 ip_rr0star *= ip_rr0star;
350 restart_flag =
false;
359 s = residual - alpha*tmp0;
366 result += alpha * p + omega * s;
367 residual = s - omega * tmp1;
370 if (residual_norm / norm_rhs_host < tag.
tolerance())
375 beta = new_ip_rr0star / ip_rr0star * alpha/omega;
376 ip_rr0star = new_ip_rr0star;
385 p = residual + beta * p;
391 tag.
error(residual_norm / norm_rhs_host);
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)
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.
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)
iterator begin()
Returns an iterator pointing to the beginning of the vector (STL like)
A tag class representing the use of no preconditioner.
Extracts the underlying context from objects.
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.
viennacl::vector< NumericT > 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.
void fast_copy(const const_vector_iterator< SCALARTYPE, ALIGNMENT > &gpu_begin, const const_vector_iterator< SCALARTYPE, ALIGNMENT > &gpu_end, CPU_ITERATOR cpu_begin)