Warning: this documentation for the development version is under construction.
00001 // Copyright (C) 2003-2009 Marc Duruflé 00002 // Copyright (C) 2001-2009 Vivien Mallet 00003 // 00004 // This file is part of the linear-algebra library Seldon, 00005 // http://seldon.sourceforge.net/. 00006 // 00007 // Seldon is free software; you can redistribute it and/or modify it under the 00008 // terms of the GNU Lesser General Public License as published by the Free 00009 // Software Foundation; either version 2.1 of the License, or (at your option) 00010 // any later version. 00011 // 00012 // Seldon is distributed in the hope that it will be useful, but WITHOUT ANY 00013 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS 00014 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for 00015 // more details. 00016 // 00017 // You should have received a copy of the GNU Lesser General Public License 00018 // along with Seldon. If not, see http://www.gnu.org/licenses/. 00019 00020 00021 #ifndef SELDON_FILE_ITERATIVE_CG_CXX 00022 00023 namespace Seldon 00024 { 00025 00027 00045 template <class Titer, class Matrix1, class Vector1, class Preconditioner> 00046 int Cg(Matrix1& A, Vector1& x, const Vector1& b, 00047 Preconditioner& M, Iteration<Titer> & iter) 00048 { 00049 const int N = A.GetM(); 00050 if (N <= 0) 00051 return 0; 00052 00053 typedef typename Vector1::value_type Complexe; 00054 Complexe rho(1), rho_1(1), alpha, beta,delta; 00055 Vector1 p(b), q(b), r(b), z(b); 00056 00057 // we initialize iter 00058 int success_init = iter.Init(b); 00059 if (success_init != 0) 00060 return iter.ErrorCode(); 00061 00062 // we compute the initial residual r = b - Ax 00063 Copy(b, r); 00064 if (!iter.IsInitGuess_Null()) 00065 MltAdd(Complexe(-1), A, x, Complexe(1), r); 00066 else 00067 x.Zero(); 00068 00069 iter.SetNumberIteration(0); 00070 // Loop until the stopping criteria are satisfied 00071 while (! iter.Finished(r)) 00072 { 00073 00074 // Preconditioning z = M^{-1} r 00075 M.Solve(A, r, z); 00076 00077 // rho = (conj(r),z) 00078 rho = DotProdConj(r, z); 00079 00080 if (rho == Complexe(0) ) 00081 { 00082 iter.Fail(1, "Cg breakdown #1"); 00083 break; 00084 } 00085 00086 if (iter.First()) 00087 Copy(z, p); 00088 else 00089 { 00090 // p = beta*p + z where beta = rho_i/rho_{i-1} 00091 beta = rho / rho_1; 00092 Mlt(beta, p); 00093 Add(Complexe(1), z, p); 00094 } 00095 00096 // matrix vector product q = A*p 00097 Mlt(A, p, q); 00098 delta = DotProdConj(p, q); 00099 if (delta == Complexe(0)) 00100 { 00101 iter.Fail(2, "Cg breakdown #2"); 00102 break; 00103 } 00104 alpha = rho / delta; 00105 00106 // x = x + alpha*p and r = r - alpha*q where alpha = rho/(bar(p),q) 00107 Add(alpha, p, x); 00108 Add(-alpha, q, r); 00109 00110 rho_1 = rho; 00111 00112 ++iter; 00113 } 00114 00115 return iter.ErrorCode(); 00116 } 00117 00118 00119 } // end namespace 00120 00121 #define SELDON_FILE_ITERATIVE_CG_CXX 00122 #endif 00123