00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021 #ifndef SELDON_FILE_ITERATIVE_MINRES_CXX
00022
00023 namespace Seldon
00024 {
00025
00027
00046 template <class Titer, class Matrix1, class Vector1, class Preconditioner>
00047 int MinRes(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 Vector1 u_old(b), u(b), r(b), v_old(b), v(b),
00056 w_old(b), w(b), z(b), w_oold(b);
00057
00058 Complexe dp, beta, ibeta, beta_old, alpha, eta, ceta;
00059 Complexe cold, coold, c, soold, sold, s, rho0, rho1, rho2, rho3;
00060
00061 int success_init = iter.Init(b);
00062 if (success_init != 0)
00063 return iter.ErrorCode();
00064
00065 Copy(b,r);
00066
00067 if (!iter.IsInitGuess_Null())
00068 MltAdd(Complexe(-1), A, x, Complexe(1), r);
00069 else
00070 x.Zero();
00071
00072 u_old.Zero(); v_old.Zero(); w_old.Zero(); w.Zero(); w_oold.Zero();
00073
00074 M.Solve(A, r, z);
00075 dp = DotProd(r, z);
00076 dp = sqrt(dp); beta = dp; eta = beta;
00077 Copy(r, v); Copy(z, u);
00078
00079 ibeta = 1.0 / beta;
00080 Mlt(ibeta, v); Mlt(ibeta, u);
00081
00082 c = 1.0; s = 0.0; cold = 1.0; sold = 0.0;
00083 Titer np = Norm2(b);
00084
00085 iter.SetNumberIteration(0);
00086
00087 while (!iter.Finished(np))
00088 {
00089
00090 Mlt(A, u, r);
00091 alpha = DotProd(r, u);
00092
00093 M.Solve(A, r, z);
00094
00095
00096
00097 Add(-alpha, v, r);
00098 Add(-alpha, u, z);
00099
00100
00101 Add(-beta, v_old, r);
00102 Add(-beta, u_old, z);
00103
00104 beta_old = beta;
00105
00106 dp = DotProd(r, z);
00107 beta = sqrt(dp);
00108
00109
00110 coold = cold; cold = c; soold = sold; sold = s;
00111
00112 rho0 = cold * alpha - coold * sold * beta_old;
00113 rho1 = sqrt(rho0*rho0 + beta*beta);
00114 rho2 = sold * alpha + coold * cold * beta_old;
00115 rho3 = soold * beta_old;
00116
00117
00118 if (rho1 == Complexe(0) )
00119 {
00120 iter.Fail(1, "Minres breakdown #1");
00121 break;
00122 }
00123 c = rho0 / rho1;
00124 s = beta / rho1;
00125
00126
00127 Copy(w_old, w_oold); Copy(w, w_old);
00128 Copy(u, w);
00129
00130 Add(-rho2, w_old, w);
00131 Add(-rho3, w_oold, w);
00132 Mlt(Complexe(1./rho1), w);
00133
00134 ceta = c*eta;
00135 Add(ceta, w, x);
00136 eta = -s*eta;
00137
00138 Copy(v, v_old); Copy(u, u_old);
00139 Copy(r, v); Copy(z, u);
00140 if (beta == Complexe(0) )
00141 {
00142 iter.Fail(2, "MinRes breakdown #2");
00143 break;
00144 }
00145 ibeta = 1.0/beta;
00146 Mlt(ibeta, v); Mlt(ibeta, u);
00147
00148
00149 np *= abs(s);
00150 ++ iter;
00151 }
00152 return iter.ErrorCode();
00153 }
00154
00155 }
00156
00157 #define SELDON_FILE_ITERATIVE_MINRES_CXX
00158 #endif