ILSSymmLqLanczos.cpp 8.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247
  1. /**
  2. * @file ILSSymmLqLanczos.cpp
  3. * @brief Iteratively solving linear equation systems with the symmetric LQ (SYMMLQ) method using Lanczos process
  4. * @author Paul Bodesheim
  5. * @date 20/01/2012 (dd/mm/yyyy)
  6. */
  7. #include <iostream>
  8. #include "ILSSymmLqLanczos.h"
  9. using namespace NICE;
  10. using namespace std;
  11. ILSSymmLqLanczos::ILSSymmLqLanczos( bool verbose, uint maxIterations, double minDelta) //, bool useFlexibleVersion )
  12. {
  13. this->minDelta = minDelta;
  14. this->maxIterations = maxIterations;
  15. this->verbose = verbose;
  16. // this->useFlexibleVersion = useFlexibleVersion;
  17. }
  18. ILSSymmLqLanczos::~ILSSymmLqLanczos()
  19. {
  20. }
  21. // void ILSSymmLqLanczos::setJacobiPreconditionerLanczos ( const Vector & jacobiPreconditioner )
  22. // {
  23. // this->jacobiPreconditioner = jacobiPreconditioner;
  24. // }
  25. int ILSSymmLqLanczos::solveLin ( const GenericMatrix & gm, const Vector & b, Vector & x )
  26. {
  27. if ( b.size() != gm.rows() ) {
  28. fthrow(Exception, "Size of vector b (" << b.size() << ") mismatches with the size of the given GenericMatrix (" << gm.rows() << ").");
  29. }
  30. if ( x.size() != gm.cols() )
  31. {
  32. x.resize(gm.cols());
  33. x.set(0.0); // bad initial solution, but whatever
  34. }
  35. // if ( verbose ) cerr << "initial solution: " << x << endl;
  36. // SYMMLQ-Method based on Lanczos vectors: implementation based on the following:
  37. //
  38. // C.C. Paige and M.A. Saunders: "Solution of sparse indefinite systems of linear equations". SIAM Journal on Numerical Analysis, p. 617--629, vol. 12, no. 4, 1975
  39. //
  40. // http://www.netlib.org/templates/templates.pdf
  41. //
  42. // declare some helpers
  43. double gamma = 0.0;
  44. double gamma_bar = 0.0;
  45. double alpha = 0.0; // alpha_j = v_j^T * A * v_j for new Lanczos vector v_j
  46. double beta = b.normL2(); // beta_1 = norm(b), in general beta_j = norm(v_j) for new Lanczos vector v_j
  47. double beta_next = 0.0; // beta_{j+1}
  48. double c_new = 0.0;
  49. double c_old = -1.0;
  50. double s_new = 0.0;
  51. double s_old = 0.0;
  52. double z_new = 0.0;
  53. double z_old = 0.0;
  54. double z_older = 0.0;
  55. double delta_new = 0.0;
  56. double epsilon_next = 0.0;
  57. // init some helping vectors
  58. Vector Av(b.size(),0.0); // Av = A * v_j
  59. Vector Ac(b.size(),0.0); // Ac = A * c_j
  60. Vector *v_new = new Vector(b.size(),0.0); // new Lanczos vector v_j
  61. Vector *v_old = 0; // Lanczos vector of the iteration before: v_{j-1}
  62. Vector *v_next = new Vector(b.size(),0.0); // Lanczos vector of the next iteration: v_{j+1}
  63. Vector *w_new = new Vector(b.size(),0.0);
  64. Vector *w_bar = new Vector(b.size(),0.0);
  65. Vector x_L (b.size(),0.0);
  66. // Vector x_C (b.size(),0.0); // x_C is a much better approximation than x_L (according to the paper mentioned above)
  67. // NOTE we store x_C in output variable x and only update this solution if the residual decreases (we are able to calculate the residual of x_C without calculating x_C)
  68. // first iteration + initialization, where b will be used as the first Lanczos vector
  69. *v_new = (1/beta)*b; // init v_1, v_1 = b / norm(b)
  70. gm.multiply(Av,*v_new); // Av = A * v_1
  71. alpha = v_new->scalarProduct(Av); // alpha_1 = v_1^T * A * v_1
  72. gamma_bar = alpha; // (gamma_bar_1 is equal to alpha_1 in ILSConjugateGradientsLanczos)
  73. *v_next = Av - (alpha*(*v_new));
  74. beta_next = v_next->normL2();
  75. v_next->normalizeL2();
  76. gamma = sqrt( (gamma_bar*gamma_bar) + (beta_next*beta_next) );
  77. c_new = gamma_bar/gamma;
  78. s_new = beta_next/gamma;
  79. z_new = beta/gamma;
  80. *w_bar = *v_new;
  81. *w_new = c_new*(*w_bar) + s_new*(*v_next);
  82. *w_bar = s_new*(*w_bar) - c_new*(*v_next);
  83. x_L = z_new*(*w_new); // first approximation of x
  84. // calculate current residual of x_C
  85. double res_x_C = (beta*beta)*(s_new*s_new)/(c_new*c_new);
  86. // store minimal residual
  87. double res_x_C_min = res_x_C;
  88. // store optimal solution x_C in output variable x instead of additional variable x_C
  89. x = x_L + (z_new/c_new)*(*w_bar); // x_C = x_L + (z_new/c_new)*(*w_bar);
  90. // calculate delta of x_L
  91. double delta_x_L = fabs(z_new) * w_new->normL2();
  92. if ( verbose ) {
  93. cerr << "ILSSymmLqLanczos: iteration 1 / " << maxIterations << endl;
  94. if ( x.size() <= 20 )
  95. cerr << "ILSSymmLqLanczos: current solution x_L: " << x_L << endl;
  96. cerr << "ILSSymmLqLanczos: delta_x_L = " << delta_x_L << endl;
  97. }
  98. // start with second iteration
  99. uint j = 2;
  100. while (j <= maxIterations )
  101. {
  102. // prepare next iteration
  103. if ( v_old == 0 ) v_old = v_new;
  104. else {
  105. delete v_old;
  106. v_old = v_new;
  107. }
  108. v_new = v_next;
  109. v_next = new Vector(b.size(),0.0);
  110. beta = beta_next;
  111. z_older = z_old;
  112. z_old = z_new;
  113. s_old = s_new;
  114. res_x_C *= (c_new*c_new);
  115. // start next iteration:
  116. // calculate next Lanczos vector v_ {j+1} based on older ones
  117. gm.multiply(Av,*v_new);
  118. alpha = v_new->scalarProduct(Av);
  119. *v_next = Av - (alpha*(*v_new)) - (beta*(*v_old)); // calculate v_{j+1}
  120. beta_next = v_next->normL2(); // calculate beta_{j+1}
  121. v_next->normalizeL2(); // normalize v_{j+1}
  122. // calculate elements of matrix L_bar_{j}
  123. gamma_bar = -c_old*s_new*beta - c_new*alpha; // calculate gamma_bar_{j}
  124. delta_new = -c_old*c_new*beta + s_new*alpha; // calculate delta_{j}
  125. //NOTE updating c_old after using it to calculate gamma_bar and delta_new is important!!
  126. c_old = c_new;
  127. // calculate helpers (equation 5.6 in the paper mentioned above)
  128. gamma = sqrt( (gamma_bar*gamma_bar) + (beta_next*beta_next) ); // calculate gamma_{j}
  129. c_new = gamma_bar/gamma; // calculate c_{j-1}
  130. s_new = beta_next/gamma; // calculate s_{j-1}
  131. // calculate next component z_{j} of vector z
  132. z_new = - (delta_new*z_old + epsilon_next*z_older)/gamma;
  133. //NOTE updating epsilon_next after using it to calculate z_new is important!!
  134. epsilon_next = s_old*beta_next; // calculate epsilon_{j+1} of matrix L_bar_{j+1}
  135. // calculate residual of current solution x_C without computing this solution x_C before
  136. res_x_C *= (s_new*s_new)/(c_new*c_new);
  137. // we only update our solution x (originally x_C ) if the residual is smaller
  138. if ( res_x_C < res_x_C_min )
  139. {
  140. x = x_L + (z_new/c_new)*(*w_bar); // x_C = x_L + (z_new/c_new)*(*w_bar); // update x
  141. res_x_C_min = res_x_C;
  142. }
  143. // calculate new vectors w_{j} and w_bar_{j+1} according to equation 5.9 of the paper mentioned above
  144. *w_new = c_new*(*w_bar) + s_new*(*v_next);
  145. *w_bar = s_new*(*w_bar) - c_new*(*v_next);
  146. x_L += z_new*(*w_new); // update x_L
  147. if ( verbose ) {
  148. cerr << "ILSSymmLqLanczos: iteration " << j << " / " << maxIterations << endl;
  149. if ( x.size() <= 20 )
  150. cerr << "ILSSymmLqLanczos: current solution x_L: " << x_L << endl;
  151. }
  152. // check convergence
  153. delta_x_L = fabs(z_new) * w_new->normL2();
  154. if ( verbose ) {
  155. cerr << "ILSSymmLqLanczos: delta_x_L = " << delta_x_L << endl;
  156. cerr << "ILSSymmLqLanczos: residual = " << res_x_C << endl;
  157. }
  158. if ( delta_x_L < minDelta ) {
  159. if ( verbose )
  160. cerr << "ILSSymmLqLanczos: small delta_x_L" << endl;
  161. break;
  162. }
  163. j++;
  164. }
  165. // if ( verbose ) {
  166. cerr << "ILSSymmLqLanczos: iterations needed: " << std::min<uint>(j,maxIterations) << endl;
  167. cerr << "ILSSymmLqLanczos: minimal residual achieved: " << res_x_C_min << endl;
  168. if ( x.size() <= 20 )
  169. cerr << "ILSSymmLqLanczos: optimal solution: " << x << endl;
  170. // }
  171. // WE DO NOT WANT TO CALCULATE THE RESIDUAL EXPLICITLY
  172. //
  173. // Vector tmp;
  174. // gm.multiply(tmp,x_C);
  175. // Vector res ( b - tmp );
  176. // double res_x_C = res.scalarProduct(res);
  177. //
  178. // gm.multiply(tmp,x_L);
  179. // res = b - tmp;
  180. // double res_x_L = res.scalarProduct(res);
  181. //
  182. // if ( res_x_L < res_x_C )
  183. // {
  184. // x = x_L;
  185. // if ( verbose )
  186. // cerr << "ILSSymmLqLanczos: x_L used with residual " << res_x_L << " < " << res_x_C << endl;
  187. //
  188. // } else
  189. // {
  190. //
  191. // x = x_C;
  192. // if ( verbose )
  193. // cerr << "ILSSymmLqLanczos: x_C used with residual " << res_x_C << " < " << res_x_L << endl;
  194. //
  195. // }
  196. delete v_new;
  197. delete v_old;
  198. delete v_next;
  199. delete w_new;
  200. delete w_bar;
  201. return 0;
  202. }