1 #ifndef VIENNACL_LINALG_CG_HPP_
2 #define VIENNACL_LINALG_CG_HPP_
64 void abs_tolerance(
double new_tol) {
if (new_tol >= 0) abs_tol_ = new_tol; }
70 unsigned int iters()
const {
return iters_taken_; }
71 void iters(
unsigned int i)
const { iters_taken_ = i; }
74 double error()
const {
return last_error_; }
76 void error(
double e)
const { last_error_ = e; }
82 unsigned int iterations_;
85 mutable unsigned int iters_taken_;
86 mutable double last_error_;
93 template<
typename VectorT,
typename PreconditionerT>
97 VectorT &
get() {
return z_; }
102 template<
typename VectorT>
105 z_handler(VectorT & residual) : presidual_(&residual){ }
106 VectorT &
get() {
return *presidual_; }
108 VectorT * presidual_;
128 template<
typename MatrixT,
typename NumericT>
134 void *monitor_data = NULL)
145 std::vector<NumericT> host_inner_prod_buffer(inner_prod_buffer.
size());
147 difference_type buffer_offset_per_vector =
static_cast<difference_type
>(buffer_size_per_vector);
154 NumericT inner_prod_rr = norm_rhs_squared;
170 inner_prod_rr = std::accumulate(host_inner_prod_buffer.begin(), host_inner_prod_buffer.begin() + buffer_offset_per_vector,
NumericT(0));
171 inner_prod_ApAp = std::accumulate(host_inner_prod_buffer.begin() + buffer_offset_per_vector, host_inner_prod_buffer.begin() + 2 * buffer_offset_per_vector,
NumericT(0));
172 inner_prod_pAp = std::accumulate(host_inner_prod_buffer.begin() + 2 * buffer_offset_per_vector, host_inner_prod_buffer.begin() + 3 * buffer_offset_per_vector,
NumericT(0));
174 if (monitor && monitor(result, std::sqrt(std::fabs(inner_prod_rr / norm_rhs_squared)), monitor_data))
179 alpha = inner_prod_rr / inner_prod_pAp;
180 beta = (alpha*alpha*inner_prod_ApAp - inner_prod_rr) / inner_prod_rr;
184 tag.
error(std::sqrt(std::fabs(inner_prod_rr) / norm_rhs_squared));
191 template<
typename NumericT>
197 void *monitor_data = NULL)
204 template<
typename NumericT>
210 void *monitor_data = NULL)
218 template<
typename NumericT>
224 void *monitor_data = NULL)
232 template<
typename NumericT>
238 void *monitor_data = NULL)
245 template<
typename NumericT>
251 void *monitor_data = NULL)
257 template<
typename MatrixT,
typename VectorT,
typename PreconditionerT>
261 PreconditionerT
const & precond,
263 void *monitor_data = NULL)
268 VectorT result = rhs;
271 VectorT residual = rhs;
274 VectorT & z = zhandler.
get();
280 CPU_NumericType alpha;
281 CPU_NumericType new_ip_rr = 0;
282 CPU_NumericType beta;
283 CPU_NumericType norm_rhs_squared = ip_rr;
284 CPU_NumericType new_ipp_rr_over_norm_rhs;
297 residual -= alpha * tmp;
301 if (static_cast<VectorT*>(&residual)==static_cast<VectorT*>(&z))
306 new_ipp_rr_over_norm_rhs = new_ip_rr / norm_rhs_squared;
307 if (monitor && monitor(result, std::sqrt(std::fabs(new_ipp_rr_over_norm_rhs)), monitor_data))
312 beta = new_ip_rr / ip_rr;
319 tag.
error(std::sqrt(std::fabs(new_ip_rr / norm_rhs_squared)));
338 template<
typename MatrixT,
typename VectorT,
typename PreconditionerT>
339 VectorT
solve(MatrixT
const & matrix, VectorT
const & rhs,
cg_tag const & tag, PreconditionerT
const & precond)
349 template<
typename IndexT,
typename NumericT,
typename PreconditionerT>
350 std::vector<NumericT>
solve(std::vector< std::map<IndexT, NumericT> >
const & A, std::vector<NumericT>
const & rhs,
cg_tag const & tag, PreconditionerT
const & precond)
360 std::vector<NumericT> result(vcl_result.
size());
371 template<
typename MatrixT,
typename VectorT>
372 VectorT
solve(MatrixT
const & matrix, VectorT
const & rhs,
cg_tag const & tag)
379 template<
typename VectorT>
385 cg_solver(
cg_tag const & tag) : tag_(tag), monitor_callback_(NULL), user_data_(NULL) {}
387 template<
typename MatrixT,
typename PreconditionerT>
388 VectorT
operator()(MatrixT
const & A, VectorT
const & b, PreconditionerT
const & precond)
const
393 mod_rhs = b - mod_rhs;
394 VectorT y =
detail::solve_impl(A, mod_rhs, tag_, precond, monitor_callback_, user_data_);
395 return init_guess_ + y;
401 template<
typename MatrixT>
402 VectorT
operator()(MatrixT
const & A, VectorT
const & b)
const
420 void set_monitor(
bool (*monitor_fun)(VectorT
const &, numeric_type,
void *),
void *user_data)
422 monitor_callback_ = monitor_fun;
423 user_data_ = user_data;
432 bool (*monitor_callback_)(VectorT
const &,
numeric_type,
void *);
Sparse matrix class using a hybrid format composed of the ELL and CSR format for storing the nonzeros...
T norm_2(std::vector< T, A > const &v1)
Generic interface for the l^2-norm. See viennacl/linalg/vector_operations.hpp for implementations...
cg_tag(double tol=1e-8, unsigned int max_iterations=300)
The constructor.
Generic size and resize functionality for different vector and matrix types.
base_type::difference_type difference_type
Generic interface for matrix-vector and matrix-matrix products. See viennacl/linalg/vector_operations...
double tolerance() const
Returns the relative tolerance.
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.
handles the no_precond case at minimal overhead
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)
void pipelined_cg_vector_update(vector_base< NumericT > &result, NumericT alpha, vector_base< NumericT > &p, vector_base< NumericT > &r, vector_base< NumericT > const &Ap, NumericT beta, vector_base< NumericT > &inner_prod_buffer)
Performs a joint vector update operation needed for an efficient pipelined CG algorithm.
unsigned int iters() const
Return the number of solver iterations:
Generic interface for the computation of inner products. See viennacl/linalg/vector_operations.hpp for implementations.
VectorT solve(MatrixT const &matrix, VectorT const &rhs, bicgstab_tag const &tag, PreconditionerT const &precond)
void pipelined_cg_prod(MatrixT const &A, vector_base< NumericT > const &p, vector_base< NumericT > &Ap, vector_base< NumericT > &inner_prod_buffer)
Performs a joint vector update operation needed for an efficient pipelined CG algorithm.
VectorT operator()(MatrixT const &A, VectorT const &b) const
Main namespace in ViennaCL. Holds all the basic types such as vector, matrix, etc. and defines operations upon them.
unsigned int max_iterations() const
Returns the maximum number of iterations.
VectorT prod(std::vector< std::vector< T, A1 >, A2 > const &matrix, VectorT const &vector)
vcl_size_t size(VectorType const &vec)
Generic routine for obtaining the size of a vector (ViennaCL, uBLAS, etc.)
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.
Implementations of incomplete factorization preconditioners. Convenience header file.
Sparse matrix class using the sliced ELLPACK with parameters C, .
z_handler(VectorT &residual)
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.
double abs_tolerance() const
Returns the absolute tolerance.
viennacl::result_of::cpu_value_type< VectorT >::type numeric_type
Generic clear functionality for different vector and matrix types.
cg_solver(cg_tag const &tag)
T::ERROR_CANNOT_DEDUCE_CPU_SCALAR_TYPE_FOR_T type
double error() const
Returns the estimated relative error at the end of the solver run.
z_handler(VectorT &residual)
void abs_tolerance(double new_tol)
Sets the absolute tolerance.
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 iters(unsigned int i) const
void error(double e) const
Sets the estimated relative error at the end of the solver run.
Helper meta function for retrieving the main RAM-based value type. Particularly important to obtain T...
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) ...
size_type size() const
Returns the length of the vector (cf. std::vector)
A tag for the conjugate gradient Used for supplying solver parameters and for dispatching the solve()...
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...
iterator end()
Returns an iterator pointing to the end of the vector (STL like)
void set_initial_guess(VectorT const &x)
Specifies an initial guess for the iterative solver.
A collection of compile time type deductions.
cg_tag const & tag() const
Returns the solver tag containing basic configuration such as tolerances, etc.
VectorT operator()(MatrixT const &A, VectorT const &b, PreconditionerT const &precond) const
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)