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_QMR_CXX 00022 00023 namespace Seldon 00024 { 00025 00027 00041 template <class Titer, class Matrix1, class Vector1, class Preconditioner> 00042 int Qmr(Matrix1& A, Vector1& x, const Vector1& b, 00043 Preconditioner& M, Iteration<Titer> & iter) 00044 { 00045 const int N = A.GetM(); 00046 if (N <= 0) 00047 return 0; 00048 00049 typedef typename Vector1::value_type Complexe; 00050 Complexe delta, ep(0), beta; 00051 Titer rho, rho_1, xi; 00052 Complexe theta_1, gamma_1; 00053 Complexe theta(0), gamma(1), eta(-1); 00054 00055 Vector1 r(b), y(b), z_tld(b); r.Zero(); 00056 Vector1 v(b), w(b), p_tld(b); 00057 Vector1 p(b), q(b), d(b), s(b); 00058 00059 // we initialize iter 00060 int success_init = iter.Init(b); 00061 if (success_init != 0) 00062 return iter.ErrorCode(); 00063 00064 // 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, v); 00072 00073 M.Solve(A, v, y); 00074 rho = Norm2(y); 00075 00076 Copy(r, w); 00077 xi = Norm2(w); 00078 00079 iter.SetNumberIteration(0); 00080 // Loop until the stopping criteria are reached 00081 while (! iter.Finished(r)) 00082 { 00083 if (rho == Titer(0)) 00084 { 00085 iter.Fail(1, "Qmr breakdown #1"); 00086 break; 00087 } 00088 if (xi == Titer(0)) 00089 { 00090 iter.Fail(2, "Qmr breakdown #2"); 00091 break; 00092 } 00093 00094 // v = v / rho 00095 // y = y / rho 00096 Mlt(Complexe(1./rho), v); 00097 Mlt(Complexe(1./rho), y); 00098 00099 // w = w / xi 00100 Mlt(Complexe(1./xi), w); 00101 00102 delta = DotProd(w, y); 00103 if (delta == Complexe(0)) 00104 { 00105 iter.Fail(3, "Qmr breakdown #3"); 00106 break; 00107 } 00108 00109 M.TransSolve(A, w, z_tld); 00110 00111 if (iter.First()) 00112 { 00113 Copy(y, p); 00114 Copy(z_tld, q); 00115 } 00116 else 00117 { 00118 // p = y - (xi delta / ep) p 00119 // q = z_tld - (rho delta / ep) q 00120 Mlt(Complexe(-(xi * delta / ep)), p); 00121 Add(Complexe(1), y, p); 00122 Mlt(Complexe(-(rho * delta / ep)), q); 00123 Add(Complexe(1), z_tld, q); 00124 } 00125 00126 // product matrix vector p_tld = A*p 00127 Mlt(A, p, p_tld); 00128 00129 ep = DotProd(q, p_tld); 00130 if (ep == Complexe(0)) 00131 { 00132 iter.Fail(4, "Qmr breakdown #4"); 00133 break; 00134 } 00135 00136 beta = ep / delta; 00137 if (beta == Complexe(0)) 00138 { 00139 iter.Fail(5, "Qmr breakdown #5"); 00140 break; 00141 } 00142 00143 // v = -beta v + p_tld 00144 Mlt(Complexe(-beta), v); 00145 Add(Complexe(1), p_tld, v); 00146 M.Solve(A, v, y); 00147 00148 rho_1 = rho; 00149 rho = Norm2(y); 00150 00151 // product matrix vector z_tld = A q 00152 Mlt(SeldonTrans, A, q, z_tld); 00153 // w = z_tld - beta*w 00154 Mlt(Complexe(-beta), w); Add(Complexe(1), z_tld, w); 00155 00156 xi = Norm2(w); 00157 00158 gamma_1 = gamma; 00159 theta_1 = theta; 00160 00161 ++iter; 00162 theta = rho / (gamma_1 * beta); 00163 gamma = Complexe(1) / sqrt(1.0 + theta * theta); 00164 00165 if (gamma == Titer(0)) 00166 { 00167 iter.Fail(6, "Qmr breakdown #6"); 00168 break; 00169 } 00170 00171 eta = -eta * rho_1 * gamma * gamma / (beta * gamma_1 * gamma_1); 00172 00173 if (iter.First()) 00174 { 00175 Copy(p, d); 00176 Mlt(eta, d); 00177 Copy(p_tld, s); 00178 Mlt(eta, s); 00179 } 00180 else 00181 { 00182 Complexe tmp = (theta_1 * theta_1 * gamma * gamma); 00183 Mlt(tmp, d); 00184 Add(eta, p, d); 00185 Mlt(tmp, s); 00186 Add(eta, p_tld, s); 00187 } 00188 Add(Complexe(1), d, x); 00189 Add(-Complexe(1), s, r); 00190 00191 ++iter; 00192 } 00193 00194 return iter.ErrorCode(); 00195 } 00196 00197 } // end namespace 00198 00199 #define SELDON_FILE_ITERATIVE_QMR_CXX 00200 #endif