ILSConjugateGradients.cpp 5.2 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207
  1. /**
  2. * @file ILSConjugateGradients.cpp
  3. * @brief abstract interface for iterative linear solvers working with GenericMatrix
  4. * @author Erik Rodner
  5. * @date 12/21/2011
  6. */
  7. #include <iostream>
  8. #include <core/basics/Timer.h>
  9. #include "ILSConjugateGradients.h"
  10. using namespace NICE;
  11. using namespace std;
  12. ILSConjugateGradients::ILSConjugateGradients( bool verbose, uint maxIterations, double minDelta, double minResidual, bool useFlexibleVersion )
  13. {
  14. this->minDelta = minDelta;
  15. this->maxIterations = maxIterations;
  16. this->verbose = verbose;
  17. this->minResidual = minResidual;
  18. this->useFlexibleVersion = useFlexibleVersion;
  19. this->timeAnalysis = false;
  20. }
  21. ILSConjugateGradients::~ILSConjugateGradients()
  22. {
  23. }
  24. void ILSConjugateGradients::setTimeAnalysis( bool timeAnalysis )
  25. {
  26. this->timeAnalysis = timeAnalysis;
  27. }
  28. void ILSConjugateGradients::setJacobiPreconditioner ( const Vector & jacobiPreconditioner )
  29. {
  30. this->jacobiPreconditioner = jacobiPreconditioner;
  31. }
  32. int ILSConjugateGradients::solveLin ( const GenericMatrix & gm, const Vector & b, Vector & x )
  33. {
  34. Timer t;
  35. if ( timeAnalysis )
  36. t.start();
  37. if ( b.size() != gm.rows() ) {
  38. fthrow(Exception, "Size of vector b (" << b.size() << ") mismatches with the size of the given GenericMatrix (" << gm.rows() << ").");
  39. }
  40. if ( x.size() != gm.cols() )
  41. {
  42. x.resize(gm.cols());
  43. x.set(0.0); // bad initial solution, but whatever
  44. }
  45. // CG-Method: http://www.netlib.org/templates/templates.pdf
  46. //
  47. Vector a;
  48. // compute r^0 = b - A*x^0
  49. gm.multiply( a, x );
  50. Vector r = b - a;
  51. if ( timeAnalysis ) {
  52. t.stop();
  53. cerr << "r = " << r << endl;
  54. cerr << "ILSConjugateGradients: TIME " << t.getSum() << " " << r.normL2() << " " << r.normInf() << endl;
  55. t.start();
  56. }
  57. Vector rold ( r.size(), 0.0 );
  58. // store optimal values
  59. double res_min = r.scalarProduct(r);
  60. Vector current_x = x;
  61. double rhoold = 0.0;
  62. Vector z ( r.size() );
  63. Vector p ( z.size() );
  64. uint i = 1;
  65. while ( i <= maxIterations )
  66. {
  67. // pre-conditioned vector, currently M=I, i.e. no pre-condition
  68. // otherwise set z = M * r
  69. if ( jacobiPreconditioner.size() != r.size() )
  70. z = r;
  71. else {
  72. // use simple Jacobi pre-conditioning
  73. for ( uint jj = 0 ; jj < z.size() ; jj++ )
  74. z[jj] = r[jj] / jacobiPreconditioner[jj];
  75. }
  76. double rho = z.scalarProduct( r );
  77. if ( verbose ) {
  78. cerr << "ILSConjugateGradients: iteration " << i << " / " << maxIterations << endl;
  79. if ( current_x.size() <= 20 )
  80. cerr << "ILSConjugateGradients: current solution " << current_x << endl;
  81. }
  82. if ( i == 1 ) {
  83. p = z;
  84. } else {
  85. double beta;
  86. if ( useFlexibleVersion ) {
  87. beta = ( rho - z.scalarProduct(rold) ) / rhoold;
  88. } else {
  89. beta = rho / rhoold;
  90. }
  91. p = z + beta * p;
  92. }
  93. Vector q ( gm.rows() );
  94. // q = A*p
  95. gm.multiply ( q, p );
  96. // sp = p^T A p
  97. // if A is next to zero this gets nasty, because we divide by sp
  98. // later on
  99. double sp = p.scalarProduct(q);
  100. if ( fabs(sp) < 1e-20 ) {
  101. // we achieved some kind of convergence, at least this
  102. // is a termination condition used in the wiki article
  103. if ( verbose )
  104. cerr << "ILSConjugateGradients: p^T*q is quite small" << endl;
  105. break;
  106. }
  107. double alpha = rho / sp;
  108. current_x = current_x + alpha * p;
  109. rold = r;
  110. r = r - alpha * q;
  111. double res = r.scalarProduct(r);
  112. double resMax = r.normInf();
  113. // store optimal x that produces smallest residual
  114. if (res < res_min) {
  115. x = current_x;
  116. res_min = res;
  117. }
  118. // check convergence
  119. double delta = fabs(alpha) * p.normL2();
  120. if ( verbose ) {
  121. cerr << "ILSConjugateGradients: delta = " << delta << " lower bound = " << minDelta << endl;
  122. cerr << "ILSConjugateGradients: residual = " << res << endl;
  123. cerr << "ILSConjugateGradients: L_inf residual = " << resMax << " lower bound = " << minResidual << endl;
  124. if (resMax < 0)
  125. {
  126. std::cerr << "WARNING: resMax is smaller than zero! " << std::endl;
  127. std::cerr << "ILSConjugateGradients: vector r: " << r << std::endl;
  128. }
  129. }
  130. if ( timeAnalysis ) {
  131. t.stop();
  132. cerr << "ILSConjugateGradients: TIME " << t.getSum() << " " << res << " " << resMax << endl;
  133. t.start();
  134. }
  135. if ( delta < minDelta ) {
  136. if ( verbose )
  137. cerr << "ILSConjugateGradients: small delta" << endl;
  138. break;
  139. }
  140. //fabs was necessary, since the IPP-implementation from normInf was not working correctly
  141. // if ( fabs(resMax) < minResidual ) {
  142. if ( resMax < minResidual ) {
  143. if ( verbose )
  144. {
  145. cerr << "ILSConjugateGradients: small residual -- resMax: "<< resMax << " minRes: " << minResidual << endl;
  146. }
  147. break;
  148. }
  149. rhoold = rho;
  150. i++;
  151. }
  152. if (verbose)
  153. {
  154. cerr << "ILSConjugateGradients: iterations needed: " << std::min<uint>(i,maxIterations) << endl;
  155. cerr << "ILSConjugateGradients: minimal residual achieved: " << res_min << endl;
  156. }
  157. if (verbose)
  158. {
  159. if ( x.size() <= 20 )
  160. cerr << "ILSConjugateGradients: optimal solution: " << x << endl;
  161. }
  162. return 0;
  163. }
  164. void ILSConjugateGradients::setVerbose(const bool& _verbose)
  165. {
  166. this->verbose = _verbose;
  167. }