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_QCGS_CXX 00022 00023 namespace Seldon 00024 { 00025 00027 00046 template <class Titer, class Matrix1, class Vector1, class Preconditioner> 00047 int QCgs(Matrix1& A, Vector1& x, const Vector1& b, 00048 Preconditioner& M, Iteration<Titer> & iter) 00049 { 00050 const int N = A.GetM(); 00051 if (N <= 0) 00052 return 0; 00053 00054 typedef typename Vector1::value_type Complexe; 00055 Complexe rho_1, rho_2, mu,nu,alpha,beta,sigma,delta; 00056 Vector1 p(b), q(b), r(b), rtilde(b), u(b), 00057 phat(b), r_qcgs(b), x_qcgs(b), v(b); 00058 00059 // we initialize iter 00060 int success_init = iter.Init(b); 00061 if (success_init != 0) 00062 return iter.ErrorCode(); 00063 00064 // we compute the residual r = b - Ax 00065 Copy(b,r); 00066 if (!iter.IsInitGuess_Null()) 00067 MltAdd(Complexe(-1), A, x, Complexe(1), r); 00068 else 00069 x.Zero(); 00070 00071 Copy(r, rtilde); 00072 rho_1 = Complexe(1); 00073 q.Zero(); p.Zero(); 00074 Copy(r, r_qcgs); Copy(x, x_qcgs); 00075 Matrix<Complexe, Symmetric, RowSymPacked> bt_b(2,2),bt_b_m1(2,2); 00076 Vector<Complexe> bt_rn(2); 00077 00078 iter.SetNumberIteration(0); 00079 // Loop until the stopping criteria are reached 00080 while (! iter.Finished(r_qcgs)) 00081 { 00082 rho_2 = DotProd(rtilde, r); 00083 00084 if (rho_1 == Complexe(0)) 00085 { 00086 iter.Fail(1, "QCgs breakdown #1"); 00087 break; 00088 } 00089 beta = rho_2/rho_1; 00090 00091 // u = r + beta*q 00092 // p = beta*(beta*p +q) + u where beta = rho_i/rho_{i-1} 00093 Copy(r, u); 00094 Add(beta, q, u); 00095 Mlt(beta, p); 00096 Add(Complexe(1), q, p); 00097 Mlt(beta, p); 00098 Add(Complexe(1), u, p); 00099 00100 // preconditioning phat = M^{-1} p 00101 M.Solve(A, p, phat); 00102 00103 // product matrix vector vhat = A*phat 00104 Mlt(A, phat, v); ++iter; 00105 sigma = DotProd(rtilde, v); 00106 if (sigma == Complexe(0)) 00107 { 00108 iter.Fail(2, "Qcgs breakdown #2"); 00109 break; 00110 } 00111 // q = u-alpha*v where alpha = rho_i/(rtilde,v) 00112 alpha = rho_2 /sigma; 00113 Copy(u, q); 00114 Add(-alpha, v, q); 00115 00116 // u = u +q 00117 Add(Complexe(1), q, u); 00118 // x = x + alpha u 00119 Add(alpha, u, x); 00120 // preconditioning phat = M^{-1} u 00121 M.Solve(A, u, phat); 00122 // product matrix vector q = A*phat 00123 Mlt(A, phat, u); 00124 00125 // r = r - alpha u 00126 Add(-alpha, u, r); 00127 // u = r_qcgs - r_n 00128 Copy(r_qcgs, u); 00129 Add(-Complexe(1), r, u); 00130 bt_b(0,0) = DotProd(u,u); 00131 bt_b(1,1) = DotProd(v,v); 00132 bt_b(1,0) = DotProd(u,v); 00133 00134 // we compute inverse of bt_b 00135 delta = bt_b(0,0)*bt_b(1,1) - bt_b(1,0)*bt_b(0,1); 00136 if (delta == Complexe(0)) 00137 { 00138 iter.Fail(3, "Qcgs breakdown #3"); 00139 break; 00140 } 00141 bt_b_m1(0,0) = bt_b(1,1)/delta; 00142 bt_b_m1(1,1) = bt_b(0,0)/delta; 00143 bt_b_m1(1,0) = -bt_b(1,0)/delta; 00144 00145 bt_rn(0) = -DotProd(u, r); bt_rn(1) = -DotProd(v, r); 00146 mu = bt_b_m1(0,0)*bt_rn(0)+bt_b_m1(0,1)*bt_rn(1); 00147 nu = bt_b_m1(1,0)*bt_rn(0)+bt_b_m1(1,1)*bt_rn(1); 00148 00149 // smoothing step 00150 // r_qcgs = r + mu (r_qcgs - r) + nu v 00151 Copy(r, r_qcgs); Add(mu, u, r_qcgs); 00152 Add(nu, v, r_qcgs); 00153 // x_qcgs = x + mu (x_qcgs - x) - nu p 00154 Mlt(mu,x_qcgs); 00155 Add(Complexe(1)-mu, x, x_qcgs); 00156 Add(-nu, p, x_qcgs); 00157 00158 rho_1 = rho_2; 00159 00160 ++iter; 00161 } 00162 00163 M.Solve(A, x_qcgs, x); 00164 return iter.ErrorCode(); 00165 } 00166 00167 } // end namespace 00168 00169 #define SELDON_FILE_ITERATIVE_QCGS_CXX 00170 #endif 00171