ViennaCL - The Vienna Computing Library
1.5.0
|
00001 #ifndef VIENNACL_LINALG_BICGSTAB_HPP_ 00002 #define VIENNACL_LINALG_BICGSTAB_HPP_ 00003 00004 /* ========================================================================= 00005 Copyright (c) 2010-2013, Institute for Microelectronics, 00006 Institute for Analysis and Scientific Computing, 00007 TU Wien. 00008 Portions of this software are copyright by UChicago Argonne, LLC. 00009 00010 ----------------- 00011 ViennaCL - The Vienna Computing Library 00012 ----------------- 00013 00014 Project Head: Karl Rupp rupp@iue.tuwien.ac.at 00015 00016 (A list of authors and contributors can be found in the PDF manual) 00017 00018 License: MIT (X11), see file LICENSE in the base directory 00019 ============================================================================= */ 00020 00025 #include <vector> 00026 #include <cmath> 00027 #include "viennacl/forwards.h" 00028 #include "viennacl/tools/tools.hpp" 00029 #include "viennacl/linalg/prod.hpp" 00030 #include "viennacl/linalg/inner_prod.hpp" 00031 #include "viennacl/linalg/norm_2.hpp" 00032 #include "viennacl/traits/clear.hpp" 00033 #include "viennacl/traits/size.hpp" 00034 #include "viennacl/meta/result_of.hpp" 00035 00036 namespace viennacl 00037 { 00038 namespace linalg 00039 { 00040 00043 class bicgstab_tag 00044 { 00045 public: 00052 bicgstab_tag(double tol = 1e-8, vcl_size_t max_iters = 400, vcl_size_t max_iters_before_restart = 200) 00053 : tol_(tol), iterations_(max_iters), iterations_before_restart_(max_iters_before_restart) {} 00054 00056 double tolerance() const { return tol_; } 00058 vcl_size_t max_iterations() const { return iterations_; } 00060 vcl_size_t max_iterations_before_restart() const { return iterations_before_restart_; } 00061 00063 vcl_size_t iters() const { return iters_taken_; } 00064 void iters(vcl_size_t i) const { iters_taken_ = i; } 00065 00067 double error() const { return last_error_; } 00069 void error(double e) const { last_error_ = e; } 00070 00071 private: 00072 double tol_; 00073 vcl_size_t iterations_; 00074 vcl_size_t iterations_before_restart_; 00075 00076 //return values from solver 00077 mutable vcl_size_t iters_taken_; 00078 mutable double last_error_; 00079 }; 00080 00081 00091 template <typename MatrixType, typename VectorType> 00092 VectorType solve(const MatrixType & matrix, VectorType const & rhs, bicgstab_tag const & tag) 00093 { 00094 typedef typename viennacl::result_of::value_type<VectorType>::type ScalarType; 00095 typedef typename viennacl::result_of::cpu_value_type<ScalarType>::type CPU_ScalarType; 00096 VectorType result = rhs; 00097 viennacl::traits::clear(result); 00098 00099 VectorType residual = rhs; 00100 VectorType p = rhs; 00101 VectorType r0star = rhs; 00102 VectorType tmp0 = rhs; 00103 VectorType tmp1 = rhs; 00104 VectorType s = rhs; 00105 00106 CPU_ScalarType norm_rhs_host = viennacl::linalg::norm_2(residual); 00107 CPU_ScalarType ip_rr0star = norm_rhs_host * norm_rhs_host; 00108 CPU_ScalarType beta; 00109 CPU_ScalarType alpha; 00110 CPU_ScalarType omega; 00111 //ScalarType inner_prod_temp; //temporary variable for inner product computation 00112 CPU_ScalarType new_ip_rr0star = 0; 00113 CPU_ScalarType residual_norm = norm_rhs_host; 00114 00115 if (norm_rhs_host == 0) //solution is zero if RHS norm is zero 00116 return result; 00117 00118 bool restart_flag = true; 00119 vcl_size_t last_restart = 0; 00120 for (vcl_size_t i = 0; i < tag.max_iterations(); ++i) 00121 { 00122 if (restart_flag) 00123 { 00124 residual = rhs; 00125 residual -= viennacl::linalg::prod(matrix, result); 00126 p = residual; 00127 r0star = residual; 00128 ip_rr0star = viennacl::linalg::norm_2(residual); 00129 ip_rr0star *= ip_rr0star; 00130 restart_flag = false; 00131 last_restart = i; 00132 } 00133 00134 tag.iters(i+1); 00135 tmp0 = viennacl::linalg::prod(matrix, p); 00136 alpha = ip_rr0star / viennacl::linalg::inner_prod(tmp0, r0star); 00137 00138 s = residual - alpha*tmp0; 00139 00140 tmp1 = viennacl::linalg::prod(matrix, s); 00141 CPU_ScalarType norm_tmp1 = viennacl::linalg::norm_2(tmp1); 00142 omega = viennacl::linalg::inner_prod(tmp1, s) / (norm_tmp1 * norm_tmp1); 00143 00144 result += alpha * p + omega * s; 00145 residual = s - omega * tmp1; 00146 00147 new_ip_rr0star = viennacl::linalg::inner_prod(residual, r0star); 00148 residual_norm = viennacl::linalg::norm_2(residual); 00149 if (std::fabs(residual_norm / norm_rhs_host) < tag.tolerance()) 00150 break; 00151 00152 beta = new_ip_rr0star / ip_rr0star * alpha/omega; 00153 ip_rr0star = new_ip_rr0star; 00154 00155 if (ip_rr0star == 0 || omega == 0 || i - last_restart > tag.max_iterations_before_restart()) //search direction degenerate. A restart might help 00156 restart_flag = true; 00157 00158 // Execution of 00159 // p = residual + beta * (p - omega*tmp0); 00160 // without introducing temporary vectors: 00161 p -= omega * tmp0; 00162 p = residual + beta * p; 00163 } 00164 00165 //store last error estimate: 00166 tag.error(residual_norm / norm_rhs_host); 00167 00168 return result; 00169 } 00170 00171 template <typename MatrixType, typename VectorType> 00172 VectorType solve(const MatrixType & matrix, VectorType const & rhs, bicgstab_tag const & tag, viennacl::linalg::no_precond) 00173 { 00174 return solve(matrix, rhs, tag); 00175 } 00176 00187 template <typename MatrixType, typename VectorType, typename PreconditionerType> 00188 VectorType solve(const MatrixType & matrix, VectorType const & rhs, bicgstab_tag const & tag, PreconditionerType const & precond) 00189 { 00190 typedef typename viennacl::result_of::value_type<VectorType>::type ScalarType; 00191 typedef typename viennacl::result_of::cpu_value_type<ScalarType>::type CPU_ScalarType; 00192 VectorType result = rhs; 00193 viennacl::traits::clear(result); 00194 00195 VectorType residual = rhs; 00196 VectorType r0star = residual; //can be chosen arbitrarily in fact 00197 VectorType tmp0 = rhs; 00198 VectorType tmp1 = rhs; 00199 VectorType s = rhs; 00200 00201 VectorType p = residual; 00202 00203 CPU_ScalarType ip_rr0star = viennacl::linalg::norm_2(residual); 00204 CPU_ScalarType norm_rhs_host = viennacl::linalg::norm_2(residual); 00205 CPU_ScalarType beta; 00206 CPU_ScalarType alpha; 00207 CPU_ScalarType omega; 00208 CPU_ScalarType new_ip_rr0star = 0; 00209 CPU_ScalarType residual_norm = norm_rhs_host; 00210 00211 if (norm_rhs_host == 0) //solution is zero if RHS norm is zero 00212 return result; 00213 00214 bool restart_flag = true; 00215 vcl_size_t last_restart = 0; 00216 for (unsigned int i = 0; i < tag.max_iterations(); ++i) 00217 { 00218 if (restart_flag) 00219 { 00220 residual = rhs; 00221 residual -= viennacl::linalg::prod(matrix, result); 00222 precond.apply(residual); 00223 p = residual; 00224 r0star = residual; 00225 ip_rr0star = viennacl::linalg::norm_2(residual); 00226 ip_rr0star *= ip_rr0star; 00227 restart_flag = false; 00228 last_restart = i; 00229 } 00230 00231 tag.iters(i+1); 00232 tmp0 = viennacl::linalg::prod(matrix, p); 00233 precond.apply(tmp0); 00234 alpha = ip_rr0star / viennacl::linalg::inner_prod(tmp0, r0star); 00235 00236 s = residual - alpha*tmp0; 00237 00238 tmp1 = viennacl::linalg::prod(matrix, s); 00239 precond.apply(tmp1); 00240 CPU_ScalarType norm_tmp1 = viennacl::linalg::norm_2(tmp1); 00241 omega = viennacl::linalg::inner_prod(tmp1, s) / (norm_tmp1 * norm_tmp1); 00242 00243 result += alpha * p + omega * s; 00244 residual = s - omega * tmp1; 00245 00246 residual_norm = viennacl::linalg::norm_2(residual); 00247 if (residual_norm / norm_rhs_host < tag.tolerance()) 00248 break; 00249 00250 new_ip_rr0star = viennacl::linalg::inner_prod(residual, r0star); 00251 00252 beta = new_ip_rr0star / ip_rr0star * alpha/omega; 00253 ip_rr0star = new_ip_rr0star; 00254 00255 if (ip_rr0star == 0 || omega == 0 || i - last_restart > tag.max_iterations_before_restart()) //search direction degenerate. A restart might help 00256 restart_flag = true; 00257 00258 // Execution of 00259 // p = residual + beta * (p - omega*tmp0); 00260 // without introducing temporary vectors: 00261 p -= omega * tmp0; 00262 p = residual + beta * p; 00263 00264 //std::cout << "Rel. Residual in current step: " << std::sqrt(std::fabs(viennacl::linalg::inner_prod(residual, residual) / norm_rhs_host)) << std::endl; 00265 } 00266 00267 //store last error estimate: 00268 tag.error(residual_norm / norm_rhs_host); 00269 00270 return result; 00271 } 00272 00273 } 00274 } 00275 00276 #endif