ILSConjugateGradientsLanczos.cpp 7.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209
  1. /**
  2. * @file ILSConjugateGradientsLanczos.cpp
  3. * @brief Iteratively solving linear equation systems with the conjugate gradients method using Lanczos process
  4. * @author Paul Bodesheim
  5. * @date 20/01/2012 (dd/mm/yyyy)
  6. */
  7. #include <iostream>
  8. #include "ILSConjugateGradientsLanczos.h"
  9. using namespace NICE;
  10. using namespace std;
  11. ILSConjugateGradientsLanczos::ILSConjugateGradientsLanczos( 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. ILSConjugateGradientsLanczos::~ILSConjugateGradientsLanczos()
  19. {
  20. }
  21. // void ILSConjugateGradientsLanczos::setJacobiPreconditionerLanczos ( const Vector & jacobiPreconditioner )
  22. // {
  23. // this->jacobiPreconditioner = jacobiPreconditioner;
  24. // }
  25. int ILSConjugateGradientsLanczos::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. // CG-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. // init some helping vectors
  43. Vector Av(b.size(),0.0); // Av = A * v_j
  44. Vector Ac(b.size(),0.0); // Ac = A * c_j
  45. Vector r(b.size(),0.0); // current residual
  46. Vector *v_new = new Vector(x.size(),0.0); // new Lanczos vector v_j
  47. Vector *v_old = new Vector(x.size(),0.0); // Lanczos vector v_{j-1} of the iteration before
  48. Vector *v_older = 0; // Lanczos vector v_{j-2} of the iteration before
  49. Vector *c_new = new Vector(x.size(),0.0); // current update vector c_j for the solution x
  50. Vector *c_old = 0; // update vector of iteration before
  51. // declare some helpers
  52. double d_new = 0; // current element of diagonal matrix D normally obtained from Cholesky factorization of tridiagonal matrix T, where T consists alpha and beta as below
  53. double d_old = 0; // corresponding element of the iteration before
  54. double l_new = 0; // current element of lower unit bidiagonal matrix L normally obtained from Cholesky factorization of tridiagonal matrix T
  55. double p_new = 0; // current element of vector p, where p is the solution of the modified linear system
  56. double p_old = 0; // corresponding element of the iteration before
  57. double alpha = 0; // alpha_j = v_j^T * A * v_j for new Lanczos vector v_j
  58. double beta = b.normL2(); // beta_1 = norm(b), in general beta_j = norm(v_j) for new Lanczos vector v_j
  59. // first iteration + initialization, where b will be used as the first Lanczos vector
  60. *v_new = (1/beta)*b; // init v_1, v_1 = b / norm(b)
  61. gm.multiply(Av,*v_new); // Av = A * v_1
  62. alpha = v_new->scalarProduct(Av); // alpha_1 = v_1^T * A * v_1
  63. d_new=alpha; // d_1 = alpha_1, d_1 is the first element of diagonal matrix D
  64. p_new = beta/d_new; // p_1 = beta_1 / d_1
  65. *c_new = *v_new; // c_1 = v_1
  66. Ac = Av; // A*c_1 = A*v_1
  67. // store current solution
  68. Vector current_x = (p_new*(*c_new)); // first approx. of x: x_1 = p_1 * c_1
  69. // calculate current residual
  70. r = b - (p_new*Ac);
  71. double res = r.scalarProduct(r);
  72. // store minimal residual
  73. double res_min = res;
  74. // store optimal solution in output variable x
  75. x = current_x;
  76. double delta_x = fabs(p_new) * c_new->normL2();
  77. if ( verbose ) {
  78. cerr << "ILSConjugateGradientsLanczos: iteration 1 / " << maxIterations << endl;
  79. if ( current_x.size() <= 20 )
  80. cerr << "ILSConjugateGradientsLanczos: current solution " << current_x << endl;
  81. cerr << "ILSConjugateGradientsLanczos: delta_x = " << delta_x << endl;
  82. cerr << "ILSConjugateGradientsLanczos: residual = " << r.scalarProduct(r) << endl;
  83. }
  84. // start with second iteration
  85. uint j = 2;
  86. while (j <= maxIterations )
  87. {
  88. // prepare d and p for next iteration
  89. d_old = d_new;
  90. p_old = p_new;
  91. // prepare vectors v_older, v_old, v_new for next iteration
  92. if ( v_older == 0) v_older = v_old;
  93. else {
  94. delete v_older;
  95. v_older = v_old;
  96. }
  97. v_old = v_new;
  98. v_new = new Vector(v_old->size(),0.0);
  99. // prepare vectors c_old, c_new for next iteration
  100. if ( c_old == 0 ) c_old = c_new;
  101. else {
  102. delete c_old;
  103. c_old = c_new;
  104. }
  105. c_new = new Vector(c_old->size(),0.0);
  106. //start next iteration:
  107. // calulate new Lanczos vector v_j based on older ones
  108. *v_new = Av - (alpha*(*v_old)) - (beta*(*v_older)); // unnormalized v_j = ( A * v_{j-1} ) - ( alpha_{j-1} * v_{j-1} ) - ( beta_{j-1} * v_{j-2} )
  109. // calculate new weight beta_j and normalize v_j
  110. beta = v_new->normL2(); // beta_j = norm(v_j)
  111. v_new->normalizeL2(); // normalize v_j
  112. // calculate new weight alpha_j
  113. gm.multiply(Av,*v_new); // Av = A * v_j
  114. alpha = v_new->scalarProduct(Av); // alpha_j = v_j^T * A * v_j
  115. // calculate l_j and d_j as the elements of the Cholesky Factorization of current tridiagonal matrix T, where T = L * D * L^T with diagonal matrix D and
  116. // lower bidiagonal matrix L; l_j and d_j are necessary for computing the new update vector c_j for the solution x and the corresponding weight
  117. l_new = beta/sqrt(d_old); // unnormalized l_j = beta_j / sqrt(d_{j-1})
  118. d_new = alpha-(l_new*l_new); // d_j = alpha_j - l_j^2
  119. l_new/=sqrt(d_old); // normalize l_j by sqrt(d_{j-1})
  120. // calculate the new weight p_j of the new update vector c_j for the solution x
  121. p_new = -p_old*l_new*d_old/d_new;
  122. // calculate the new update vector c_j for the solution x
  123. *c_new = *v_new - (l_new*(*c_old));
  124. // calculate new residual vector
  125. Ac = Av - (l_new*Ac);
  126. r-=p_new*Ac;
  127. res = r.scalarProduct(r);
  128. // update solution x
  129. current_x+=(p_new*(*c_new));
  130. if ( verbose ) {
  131. cerr << "ILSConjugateGradientsLanczos: iteration " << j << " / " << maxIterations << endl;
  132. if ( current_x.size() <= 20 )
  133. cerr << "ILSConjugateGradientsLanczos: current solution " << current_x << endl;
  134. }
  135. // store optimal x that produces smallest residual
  136. if (res < res_min) {
  137. x = current_x;
  138. res_min = res;
  139. }
  140. // check convergence
  141. delta_x = fabs(p_new) * c_new->normL2();
  142. if ( verbose ) {
  143. cerr << "ILSConjugateGradientsLanczos: delta_x = " << delta_x << endl;
  144. cerr << "ILSConjugateGradientsLanczos: residual = " << r.scalarProduct(r) << endl;
  145. }
  146. if ( delta_x < minDelta ) {
  147. if ( verbose )
  148. cerr << "ILSConjugateGradientsLanczos: small delta_x" << endl;
  149. break;
  150. }
  151. j++;
  152. }
  153. // if ( verbose ) {
  154. cerr << "ILSConjugateGradientsLanczos: iterations needed: " << std::min<uint>(j,maxIterations) << endl;
  155. cerr << "ILSConjugateGradientsLanczos: minimal residual achieved: " << res_min << endl;
  156. if ( x.size() <= 20 )
  157. cerr << "ILSConjugateGradientsLanczos: optimal solution: " << x << endl;
  158. // }
  159. delete v_new;
  160. delete v_old;
  161. delete v_older;
  162. delete c_new;
  163. delete c_old;
  164. return 0;
  165. }