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_CGS_CXX 00022 00023 namespace Seldon 00024 { 00025 00027 00045 template <class Titer, class Matrix1, class Vector1, class Preconditioner> 00046 int Cgs(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_2(0), alpha, beta, delta; 00055 Vector1 p(b), phat(b), q(b), qhat(b), vhat(b), u(b), uhat(b), 00056 r(b), rtilde(b); 00057 00058 // we initialize iter 00059 int success_init = iter.Init(b); 00060 if (success_init != 0) 00061 return iter.ErrorCode(); 00062 00063 // we compute the initial residual r = b - Ax 00064 Copy(b,r); 00065 if (!iter.IsInitGuess_Null()) 00066 MltAdd(Complexe(-1), A, x, Complexe(1), r); 00067 else 00068 x.Zero(); 00069 00070 Copy(r, rtilde); 00071 00072 iter.SetNumberIteration(0); 00073 // Loop until the stopping criteria are reached 00074 while (! iter.Finished(r)) 00075 { 00076 rho_1 = DotProd(rtilde, r); 00077 00078 if (rho_1 == Complexe(0)) 00079 { 00080 iter.Fail(1, "Cgs breakdown #1"); 00081 break; 00082 } 00083 00084 if (iter.First()) 00085 { 00086 Copy(r, u); 00087 Copy(u, p); 00088 } 00089 else 00090 { 00091 // u = r + beta*q 00092 // p = beta*(beta*p +q) + u where beta = rho_i/rho_{i-1} 00093 beta = rho_1 / rho_2; 00094 Copy(r, u); 00095 Add(beta, q, u); 00096 Mlt(beta, p); 00097 Add(Complexe(1), q, p); 00098 Mlt(beta, p); 00099 Add(Complexe(1), u, p); 00100 } 00101 00102 // preconditioning phat = M^{-1} p 00103 M.Solve(A, p, phat); 00104 00105 // matrix vector product vhat = A*phat 00106 Mlt(A, phat, vhat); ++iter; 00107 delta = DotProd(rtilde, vhat); 00108 if (delta == Complexe(0)) 00109 { 00110 iter.Fail(2, "Cgs breakdown #2"); 00111 break; 00112 } 00113 // q = u-alpha*vhat where alpha = rho_i/(rtilde,vhat) 00114 alpha = rho_1 /delta; 00115 Copy(u,q); 00116 Add(-alpha, vhat, q); 00117 00118 // u =u+q 00119 Add(Complexe(1), q, u); 00120 M.Solve(A, u, uhat); 00121 00122 Add(alpha, uhat, x); 00123 Mlt(A, uhat, qhat); 00124 Add(-alpha, qhat, r); 00125 00126 rho_2 = rho_1; 00127 00128 ++iter; 00129 } 00130 00131 return iter.ErrorCode(); 00132 } 00133 00134 } // end namespace 00135 00136 #define SELDON_FILE_ITERATIVE_CGS_CXX 00137 #endif 00138