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_LSQR_CXX 00021 00022 namespace Seldon 00023 { 00024 00026 00040 template <class Titer, class Matrix1, class Vector1, class Preconditioner> 00041 int Lsqr(Matrix1& A, Vector1& x, const Vector1& b, 00042 Preconditioner& M, Iteration<Titer> & iter) 00043 { 00044 const int N = A.GetM(); 00045 if (N <= 0) 00046 return 0; 00047 00048 typedef typename Vector1::value_type Complexe; 00049 Complexe rho, rho_bar, phi, phi_bar, theta, c, s, tmp; 00050 Titer beta, alpha, rnorm; 00051 Vector1 v(b), v1(b), u(b), u1(b), w(b); 00052 00053 int success_init = iter.Init(b); 00054 if (success_init != 0) 00055 return iter.ErrorCode(); 00056 00057 Copy(b, u); 00058 if (!iter.IsInitGuess_Null()) 00059 MltAdd(Complexe(-1), A, x, Complexe(1), u); 00060 else 00061 x.Zero(); 00062 00063 rnorm = Norm2(u); 00064 00065 Copy(b, u); 00066 beta = Norm2(u); 00067 tmp = 1.0/beta; Mlt(tmp, u); 00068 // matrix vector product 00069 Mlt(SeldonTrans,A, u, v); 00070 alpha = Norm2(v); 00071 tmp = 1.0/alpha; Mlt(tmp, v); 00072 00073 Copy(v,w); x.Zero(); 00074 00075 phi_bar = beta; rho_bar = alpha; 00076 00077 iter.SetNumberIteration(0); 00078 // Loop until the stopping criteria are satisfied 00079 while (! iter.Finished(rnorm)) 00080 { 00081 // matrix vector product u1 = A*v 00082 Mlt(A, v, u1); 00083 // u1 = u1 - alpha*u 00084 Add(-alpha, u, u1); 00085 beta = Norm2(u1); 00086 if (beta == Complexe(0) ) 00087 { 00088 iter.Fail(1, "Lsqr breakdown #1"); 00089 break; 00090 } 00091 tmp = 1.0/beta; Mlt(tmp, u1); 00092 00093 // matrix vector product v1 = A^t u1 00094 Mlt(SeldonTrans, A, u1, v1); 00095 // v1 = v1 - beta*v 00096 Add(-beta, v, v1); 00097 alpha = Norm2(v1); 00098 if (alpha == Complexe(0) ) 00099 { 00100 iter.Fail(2, "Lsqr breakdown #2"); 00101 break; 00102 } 00103 tmp = 1.0/alpha; Mlt(tmp, v1); 00104 00105 rho = sqrt(rho_bar*rho_bar+beta*beta); 00106 if (rho == Complexe(0) ) 00107 { 00108 iter.Fail(3, "Lsqr breakdown #3"); 00109 break; 00110 } 00111 00112 c = rho_bar/rho; 00113 s = beta / rho; 00114 theta = s*alpha; 00115 rho_bar = -c * alpha; 00116 phi = c * phi_bar; 00117 phi_bar = s * phi_bar; 00118 00119 // x = x + (phi/rho) w 00120 tmp = phi/rho; 00121 Add(tmp, w, x); 00122 // w = v1 - (theta/rho) w 00123 tmp = -theta/rho; 00124 Mlt(tmp,w); Add(Complexe(1), v1, w); 00125 00126 rnorm = abs(phi_bar); 00127 00128 Swap(u1,u); Swap(v1,v); 00129 00130 ++iter; 00131 ++iter; 00132 } 00133 00134 return iter.ErrorCode(); 00135 } 00136 00137 } // end namespace 00138 00139 #define SELDON_FILE_ITERATIVE_LSQR_CXX 00140 #endif