Warning: this documentation for the development version is under construction.
00001 // Copyright (C) 2003-2009 Marc Duruflé 00002 // 00003 // This file is part of the linear-algebra library Seldon, 00004 // http://seldon.sourceforge.net/. 00005 // 00006 // Seldon is free software; you can redistribute it and/or modify it under the 00007 // terms of the GNU Lesser General Public License as published by the Free 00008 // Software Foundation; either version 2.1 of the License, or (at your option) 00009 // any later version. 00010 // 00011 // Seldon is distributed in the hope that it will be useful, but WITHOUT ANY 00012 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS 00013 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for 00014 // more details. 00015 // 00016 // You should have received a copy of the GNU Lesser General Public License 00017 // along with Seldon. If not, see http://www.gnu.org/licenses/. 00018 00019 00020 #ifndef SELDON_FILE_ITERATIVE_SYMMLQ_CXX 00021 00022 namespace Seldon 00023 { 00024 00026 00045 template <class Titer, class Matrix1, class Vector1, class Preconditioner> 00046 int Symmlq(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 alpha, beta, ibeta, beta_old, beta1, 00055 ceta(0), ceta_oold, ceta_old, ceta_bar; 00056 Complexe c, cold, s, sold, coold, soold, rho0, rho1, rho2, rho3, dp; 00057 00058 Vector1 r(b), z(b), u(b), v(b), w(b), u_old(b), v_old(b), w_bar(b); 00059 00060 Titer np, s_prod; 00061 u_old.Zero(); v_old.Zero(); w.Zero(); w_bar.Zero(); 00062 00063 int success_init = iter.Init(b); 00064 if (success_init != 0) 00065 return iter.ErrorCode(); 00066 00067 Copy(b, r); 00068 // r = b - A x 00069 if (!iter.IsInitGuess_Null()) 00070 MltAdd(Complexe(-1), A, x, Complexe(1), r); 00071 else 00072 x.Zero(); 00073 00074 ceta_oold = 0.0; ceta_old = 0.0; 00075 c = 1.0; cold = 1.0; s = 0.0; sold = 0.0; 00076 M.Solve(A, r, z); 00077 00078 dp = DotProd(r, z); 00079 dp = sqrt(dp); beta = dp; beta1 = beta; 00080 s_prod = abs(beta1); 00081 00082 Copy(r, v); Copy(z, u); 00083 ibeta = 1.0/beta; 00084 Mlt(ibeta, v); Mlt(ibeta, u); 00085 Copy(u, w_bar); 00086 np = Norm2(b); 00087 00088 iter.SetNumberIteration(0); 00089 // Loop until the stopping criteria are satisfied 00090 while (!iter.Finished(np)) 00091 { 00092 // update 00093 if (!iter.First()) 00094 { 00095 Copy(v, v_old); Copy(u, u_old); 00096 ibeta = 1.0/beta; 00097 // v = ibeta r 00098 // u = ibeta z 00099 Copy(r, v); Mlt(ibeta, v); 00100 Copy(z, u); Mlt(ibeta, u); 00101 // w = c*w_bar + s*u 00102 Copy(w_bar, w); Mlt(c, w); Add(s, u, w); 00103 // w_bar = -s*w_bar + c*u 00104 Mlt(Complexe(-s),w_bar); Add(c,u,w_bar); 00105 // x = x+ceta*w 00106 Add(ceta,w,x); 00107 00108 ceta_oold = ceta_old; 00109 ceta_old = ceta; 00110 } 00111 00112 // product matrix vector r = A u 00113 Mlt(A, u, r); 00114 alpha = DotProd(u,r); 00115 // preconditioning 00116 M.Solve(A, r,z); 00117 00118 // r = r - alpha*v 00119 // z = z - alpha*u 00120 Add(-alpha,v,r); 00121 Add(-alpha,u,z); 00122 00123 // r = r - beta*v_old 00124 // z = z - beta*u_old 00125 Add(-beta,v_old,r); 00126 Add(-beta,u_old,z); 00127 00128 beta_old = beta; 00129 dp = DotProd(r,z); 00130 beta = sqrt(dp); 00131 00132 // QR factorization 00133 coold = cold; cold = c; soold = sold; sold = s; 00134 rho0 = cold * alpha - coold * sold * beta_old; // gamma_bar 00135 rho1 = sqrt(rho0*rho0 + beta*beta); // gamma 00136 rho2 = sold * alpha + coold * cold * beta_old; // delta 00137 rho3 = soold * beta_old; // epsilon 00138 00139 // Givens rotation 00140 c = rho0 / rho1; s = beta / rho1; 00141 00142 if (iter.First()) 00143 ceta = beta1/rho1; 00144 else 00145 ceta = -(rho2*ceta_old + rho3*ceta_oold)/rho1; 00146 00147 s_prod *= abs(s); 00148 if (c == Complexe(0)) 00149 np = s_prod*1e16; 00150 else 00151 np = s_prod/abs(c); 00152 00153 ++iter; 00154 } 00155 if (c == Complexe(0)) 00156 ceta_bar = ceta*1e15; 00157 else 00158 ceta_bar = ceta/c; 00159 00160 Add(ceta_bar,w_bar,x); 00161 00162 return iter.ErrorCode(); 00163 } 00164 00165 } // end namespace 00166 00167 #define SELDON_FILE_ITERATIVE_SYMMLQ_CXX 00168 #endif