]> gitweb.michael.orlitzky.com - octave.git/blobdiff - optimization/preconditioned_conjugate_gradient_method.m
Add a second reference for the PCGM and make it more resistant to accumulated roundof...
[octave.git] / optimization / preconditioned_conjugate_gradient_method.m
index 63943482c8c6dd1b47d916d9bdaf400521d6b992..4e68ccb583e7fab2bc8cf9f34ce170eeab7fe3a2 100644 (file)
@@ -1,8 +1,8 @@
-function [x, k] = preconditioned_conjugate_gradient_method(Q,
-                                                          M,
-                                                          b,
-                                                          x0,
-                                                          tolerance,
+function [x, k] = preconditioned_conjugate_gradient_method(Q, ...
+                                                          M, ...
+                                                          b, ...
+                                                          x0, ...
+                                                          tolerance, ...
                                                           max_iterations)
   %
   % Solve,
@@ -39,7 +39,7 @@ function [x, k] = preconditioned_conjugate_gradient_method(Q,
   %
   % OUTPUT:
   %
-  %   - ``x`` - The solution to Qx=b.
+  %   - ``x`` - The computed solution to Qx=b.
   %
   %   - ``k`` - The ending value of k; that is, the number of
   %   iterations that were performed.
@@ -52,14 +52,24 @@ function [x, k] = preconditioned_conjugate_gradient_method(Q,
   % Conjugate-Gradient Method", we are supposed to define
   % d_{0} = -z_{0}, not -r_{0} as written.
   %
+  % The rather verbose name of this function was chosen to avoid
+  % conflicts with other implementations.
+  %
   % REFERENCES:
   %
   %   1. Guler, Osman. Foundations of Optimization. New York, Springer,
-  %   2010.
+  %      2010.
+  %
+  %   2. Shewchuk, Jonathan Richard. An Introduction to the Conjugate
+  %      Gradient Method Without the Agonizing Pain, Edition 1.25.
+  %      August 4, 1994.
   %
 
+  % We use this in the inner loop.
+  sqrt_n = floor(sqrt(length(x0)));
+
   % Set k=0 first, that way the references to xk,rk,zk,dk which
-  % immediately follow correspond to x0,r0,z0,d0 respectively.
+  % immediately follow correspond (semantically) to x0,r0,z0,d0.
   k = 0;
 
   xk = x0;
@@ -68,9 +78,11 @@ function [x, k] = preconditioned_conjugate_gradient_method(Q,
   dk = -zk;
 
   for k = [ 0 : max_iterations ]
+
     if (norm(rk) < tolerance)
-       x = xk;
-       return;
+      % Check our stopping condition. This should catch the k=0 case.
+      x = xk;
+      return;
     end
 
     % Used twice, avoid recomputation.
@@ -82,7 +94,17 @@ function [x, k] = preconditioned_conjugate_gradient_method(Q,
 
     alpha_k = rkzk/(dk' * Qdk);
     x_next = xk + (alpha_k * dk);
-    r_next = rk + (alpha_k * Qdk);
+
+    % The recursive definition of r_next is prone to accumulate
+    % roundoff error. When sqrt(n) divides k, we recompute the
+    % residual to minimize this error. This modification is due to the
+    % second reference.
+    if (mod(k, sqrt_n) == 0)
+      r_next = Q*x_next - b;
+    else
+      r_next = rk + (alpha_k * Qdk);
+    end
+
     z_next = M \ r_next;
     beta_next = (r_next' * z_next)/rkzk;
     d_next = -z_next + beta_next*dk;
@@ -94,5 +116,7 @@ function [x, k] = preconditioned_conjugate_gradient_method(Q,
     dk = d_next;
   end
 
+  % The algorithm didn't converge, but we still want to return the
+  % terminal value of xk.
   x = xk;
 end