View Javadoc

1   /*
2    * Licensed to the Apache Software Foundation (ASF) under one or more
3    * contributor license agreements.  See the NOTICE file distributed with
4    * this work for additional information regarding copyright ownership.
5    * The ASF licenses this file to You under the Apache License, Version 2.0
6    * (the "License"); you may not use this file except in compliance with
7    * the License.  You may obtain a copy of the License at
8    *
9    *      http://www.apache.org/licenses/LICENSE-2.0
10   *
11   * Unless required by applicable law or agreed to in writing, software
12   * distributed under the License is distributed on an "AS IS" BASIS,
13   * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14   * See the License for the specific language governing permissions and
15   * limitations under the License.
16   */
17  package org.apache.commons.math.estimation;
18  
19  import java.io.Serializable;
20  import java.util.Arrays;
21  
22  
23  /** 
24   * This class solves a least squares problem.
25   *
26   * <p>This implementation <em>should</em> work even for over-determined systems
27   * (i.e. systems having more variables than equations). Over-determined systems
28   * are solved by ignoring the variables which have the smallest impact according
29   * to their jacobian column norm. Only the rank of the matrix and some loop bounds
30   * are changed to implement this. This feature has undergone only basic testing
31   * for now and should still be considered experimental.</p>
32   *
33   * <p>The resolution engine is a simple translation of the MINPACK <a
34   * href="http://www.netlib.org/minpack/lmder.f">lmder</a> routine with minor
35   * changes. The changes include the over-determined resolution and the Q.R.
36   * decomposition which has been rewritten following the algorithm described in the
37   * P. Lascaux and R. Theodor book <i>Analyse num&eacute;rique matricielle
38   * appliqu&eacute;e &agrave; l'art de l'ing&eacute;nieur</i>, Masson 1986. The
39   * redistribution policy for MINPACK is available <a
40   * href="http://www.netlib.org/minpack/disclaimer">here</a>, for convenience, it
41   * is reproduced below.</p>
42   *
43   * <table border="0" width="80%" cellpadding="10" align="center" bgcolor="#E0E0E0">
44   * <tr><td>
45   *    Minpack Copyright Notice (1999) University of Chicago.
46   *    All rights reserved
47   * </td></tr>
48   * <tr><td>
49   * Redistribution and use in source and binary forms, with or without
50   * modification, are permitted provided that the following conditions
51   * are met:
52   * <ol>
53   *  <li>Redistributions of source code must retain the above copyright
54   *      notice, this list of conditions and the following disclaimer.</li>
55   * <li>Redistributions in binary form must reproduce the above
56   *     copyright notice, this list of conditions and the following
57   *     disclaimer in the documentation and/or other materials provided
58   *     with the distribution.</li>
59   * <li>The end-user documentation included with the redistribution, if any,
60   *     must include the following acknowledgment:
61   *     <code>This product includes software developed by the University of
62   *           Chicago, as Operator of Argonne National Laboratory.</code>
63   *     Alternately, this acknowledgment may appear in the software itself,
64   *     if and wherever such third-party acknowledgments normally appear.</li>
65   * <li><strong>WARRANTY DISCLAIMER. THE SOFTWARE IS SUPPLIED "AS IS"
66   *     WITHOUT WARRANTY OF ANY KIND. THE COPYRIGHT HOLDER, THE
67   *     UNITED STATES, THE UNITED STATES DEPARTMENT OF ENERGY, AND
68   *     THEIR EMPLOYEES: (1) DISCLAIM ANY WARRANTIES, EXPRESS OR
69   *     IMPLIED, INCLUDING BUT NOT LIMITED TO ANY IMPLIED WARRANTIES
70   *     OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE, TITLE
71   *     OR NON-INFRINGEMENT, (2) DO NOT ASSUME ANY LEGAL LIABILITY
72   *     OR RESPONSIBILITY FOR THE ACCURACY, COMPLETENESS, OR
73   *     USEFULNESS OF THE SOFTWARE, (3) DO NOT REPRESENT THAT USE OF
74   *     THE SOFTWARE WOULD NOT INFRINGE PRIVATELY OWNED RIGHTS, (4)
75   *     DO NOT WARRANT THAT THE SOFTWARE WILL FUNCTION
76   *     UNINTERRUPTED, THAT IT IS ERROR-FREE OR THAT ANY ERRORS WILL
77   *     BE CORRECTED.</strong></li>
78   * <li><strong>LIMITATION OF LIABILITY. IN NO EVENT WILL THE COPYRIGHT
79   *     HOLDER, THE UNITED STATES, THE UNITED STATES DEPARTMENT OF
80   *     ENERGY, OR THEIR EMPLOYEES: BE LIABLE FOR ANY INDIRECT,
81   *     INCIDENTAL, CONSEQUENTIAL, SPECIAL OR PUNITIVE DAMAGES OF
82   *     ANY KIND OR NATURE, INCLUDING BUT NOT LIMITED TO LOSS OF
83   *     PROFITS OR LOSS OF DATA, FOR ANY REASON WHATSOEVER, WHETHER
84   *     SUCH LIABILITY IS ASSERTED ON THE BASIS OF CONTRACT, TORT
85   *     (INCLUDING NEGLIGENCE OR STRICT LIABILITY), OR OTHERWISE,
86   *     EVEN IF ANY OF SAID PARTIES HAS BEEN WARNED OF THE
87   *     POSSIBILITY OF SUCH LOSS OR DAMAGES.</strong></li>
88   * <ol></td></tr>
89   * </table>
90  
91   * @author Argonne National Laboratory. MINPACK project. March 1980 (original fortran)
92   * @author Burton S. Garbow (original fortran)
93   * @author Kenneth E. Hillstrom (original fortran)
94   * @author Jorge J. More (original fortran)
95  
96   * @version $Revision: 620312 $ $Date: 2008-02-10 12:28:59 -0700 (Sun, 10 Feb 2008) $
97   * @since 1.2
98   *
99   */
100 public class LevenbergMarquardtEstimator extends AbstractEstimator implements Serializable {
101 
102   /** 
103    * Build an estimator for least squares problems.
104    * <p>The default values for the algorithm settings are:
105    *   <ul>
106    *    <li>{@link #setInitialStepBoundFactor initial step bound factor}: 100.0</li>
107    *    <li>{@link #setMaxCostEval maximal cost evaluations}: 1000</li>
108    *    <li>{@link #setCostRelativeTolerance cost relative tolerance}: 1.0e-10</li>
109    *    <li>{@link #setParRelativeTolerance parameters relative tolerance}: 1.0e-10</li>
110    *    <li>{@link #setOrthoTolerance orthogonality tolerance}: 1.0e-10</li>
111    *   </ul>
112    * </p>
113    */
114   public LevenbergMarquardtEstimator() {
115 
116     // set up the superclass with a default  max cost evaluations setting
117     setMaxCostEval(1000);
118 
119     // default values for the tuning parameters
120     setInitialStepBoundFactor(100.0);
121     setCostRelativeTolerance(1.0e-10);
122     setParRelativeTolerance(1.0e-10);
123     setOrthoTolerance(1.0e-10);
124 
125   }
126 
127   /** 
128    * Set the positive input variable used in determining the initial step bound.
129    * This bound is set to the product of initialStepBoundFactor and the euclidean norm of diag*x if nonzero,
130    * or else to initialStepBoundFactor itself. In most cases factor should lie
131    * in the interval (0.1, 100.0). 100.0 is a generally recommended value
132    * 
133    * @param initialStepBoundFactor initial step bound factor
134    * @see #estimate
135    */
136   public void setInitialStepBoundFactor(double initialStepBoundFactor) {
137     this.initialStepBoundFactor = initialStepBoundFactor;
138   }
139 
140   /** 
141    * Set the desired relative error in the sum of squares.
142    * 
143    * @param costRelativeTolerance desired relative error in the sum of squares
144    * @see #estimate
145    */
146   public void setCostRelativeTolerance(double costRelativeTolerance) {
147     this.costRelativeTolerance = costRelativeTolerance;
148   }
149 
150   /** 
151    * Set the desired relative error in the approximate solution parameters.
152    * 
153    * @param parRelativeTolerance desired relative error
154    * in the approximate solution parameters
155    * @see #estimate
156    */
157   public void setParRelativeTolerance(double parRelativeTolerance) {
158     this.parRelativeTolerance = parRelativeTolerance;
159   }
160 
161   /** 
162    * Set the desired max cosine on the orthogonality.
163    * 
164    * @param orthoTolerance desired max cosine on the orthogonality
165    * between the function vector and the columns of the jacobian
166    * @see #estimate
167    */
168   public void setOrthoTolerance(double orthoTolerance) {
169     this.orthoTolerance = orthoTolerance;
170   }
171 
172   /** 
173    * Solve an estimation problem using the Levenberg-Marquardt algorithm.
174    * <p>The algorithm used is a modified Levenberg-Marquardt one, based
175    * on the MINPACK <a href="http://www.netlib.org/minpack/lmder.f">lmder</a>
176    * routine. The algorithm settings must have been set up before this method
177    * is called with the {@link #setInitialStepBoundFactor},
178    * {@link #setMaxCostEval}, {@link #setCostRelativeTolerance},
179    * {@link #setParRelativeTolerance} and {@link #setOrthoTolerance} methods.
180    * If these methods have not been called, the default values set up by the
181    * {@link #LevenbergMarquardtEstimator() constructor} will be used.</p>
182    * <p>The authors of the original fortran function are:</p>
183    * <ul>
184    *   <li>Argonne National Laboratory. MINPACK project. March 1980</li>
185    *   <li>Burton  S. Garbow</li>
186    *   <li>Kenneth E. Hillstrom</li>
187    *   <li>Jorge   J. More</li>
188    *   </ul>
189    * <p>Luc Maisonobe did the Java translation.</p>
190    * 
191    * @param problem estimation problem to solve
192    * @exception EstimationException if convergence cannot be
193    * reached with the specified algorithm settings or if there are more variables
194    * than equations
195    * @see #setInitialStepBoundFactor
196    * @see #setCostRelativeTolerance
197    * @see #setParRelativeTolerance
198    * @see #setOrthoTolerance
199    */
200   public void estimate(EstimationProblem problem)
201     throws EstimationException {
202 
203     initializeEstimate(problem);
204 
205     // arrays shared with the other private methods
206     solvedCols  = Math.min(rows, cols);
207     diagR       = new double[cols];
208     jacNorm     = new double[cols];
209     beta        = new double[cols];
210     permutation = new int[cols];
211     lmDir       = new double[cols];
212 
213     // local variables
214     double   delta   = 0, xNorm = 0;
215     double[] diag    = new double[cols];
216     double[] oldX    = new double[cols];
217     double[] oldRes  = new double[rows];
218     double[] work1   = new double[cols];
219     double[] work2   = new double[cols];
220     double[] work3   = new double[cols];
221 
222     // evaluate the function at the starting point and calculate its norm
223     updateResidualsAndCost();
224     
225     // outer loop
226     lmPar = 0;
227     boolean firstIteration = true;
228     while (true) {
229 
230       // compute the Q.R. decomposition of the jacobian matrix
231       updateJacobian();
232       qrDecomposition();
233 
234       // compute Qt.res
235       qTy(residuals);
236 
237       // now we don't need Q anymore,
238       // so let jacobian contain the R matrix with its diagonal elements
239       for (int k = 0; k < solvedCols; ++k) {
240         int pk = permutation[k];
241         jacobian[k * cols + pk] = diagR[pk];
242       }
243 
244       if (firstIteration) {
245 
246         // scale the variables according to the norms of the columns
247         // of the initial jacobian
248         xNorm = 0;
249         for (int k = 0; k < cols; ++k) {
250           double dk = jacNorm[k];
251           if (dk == 0) {
252             dk = 1.0;
253           }
254           double xk = dk * parameters[k].getEstimate();
255           xNorm  += xk * xk;
256           diag[k] = dk;
257         }
258         xNorm = Math.sqrt(xNorm);
259         
260         // initialize the step bound delta
261         delta = (xNorm == 0) ? initialStepBoundFactor : (initialStepBoundFactor * xNorm);
262  
263       }
264 
265       // check orthogonality between function vector and jacobian columns
266       double maxCosine = 0;
267       if (cost != 0) {
268         for (int j = 0; j < solvedCols; ++j) {
269           int    pj = permutation[j];
270           double s  = jacNorm[pj];
271           if (s != 0) {
272             double sum = 0;
273             for (int i = 0, index = pj; i <= j; ++i, index += cols) {
274               sum += jacobian[index] * residuals[i];
275             }
276             maxCosine = Math.max(maxCosine, Math.abs(sum) / (s * cost));
277           }
278         }
279       }
280       if (maxCosine <= orthoTolerance) {
281         return;
282       }
283 
284       // rescale if necessary
285       for (int j = 0; j < cols; ++j) {
286         diag[j] = Math.max(diag[j], jacNorm[j]);
287       }
288 
289       // inner loop
290       for (double ratio = 0; ratio < 1.0e-4;) {
291 
292         // save the state
293         for (int j = 0; j < solvedCols; ++j) {
294           int pj = permutation[j];
295           oldX[pj] = parameters[pj].getEstimate();
296         }
297         double previousCost = cost;
298         double[] tmpVec = residuals;
299         residuals = oldRes;
300         oldRes    = tmpVec;
301         
302         // determine the Levenberg-Marquardt parameter
303         determineLMParameter(oldRes, delta, diag, work1, work2, work3);
304 
305         // compute the new point and the norm of the evolution direction
306         double lmNorm = 0;
307         for (int j = 0; j < solvedCols; ++j) {
308           int pj = permutation[j];
309           lmDir[pj] = -lmDir[pj];
310           parameters[pj].setEstimate(oldX[pj] + lmDir[pj]);
311           double s = diag[pj] * lmDir[pj];
312           lmNorm  += s * s;
313         }
314         lmNorm = Math.sqrt(lmNorm);
315 
316         // on the first iteration, adjust the initial step bound.
317         if (firstIteration) {
318           delta = Math.min(delta, lmNorm);
319         }
320 
321         // evaluate the function at x + p and calculate its norm
322         updateResidualsAndCost();
323 
324         // compute the scaled actual reduction
325         double actRed = -1.0;
326         if (0.1 * cost < previousCost) {
327           double r = cost / previousCost;
328           actRed = 1.0 - r * r;
329         }
330 
331         // compute the scaled predicted reduction
332         // and the scaled directional derivative
333         for (int j = 0; j < solvedCols; ++j) {
334           int pj = permutation[j];
335           double dirJ = lmDir[pj];
336           work1[j] = 0;
337           for (int i = 0, index = pj; i <= j; ++i, index += cols) {
338             work1[i] += jacobian[index] * dirJ;
339           }
340         }
341         double coeff1 = 0;
342         for (int j = 0; j < solvedCols; ++j) {
343          coeff1 += work1[j] * work1[j];
344         }
345         double pc2 = previousCost * previousCost;
346         coeff1 = coeff1 / pc2;
347         double coeff2 = lmPar * lmNorm * lmNorm / pc2;
348         double preRed = coeff1 + 2 * coeff2;
349         double dirDer = -(coeff1 + coeff2);
350 
351         // ratio of the actual to the predicted reduction
352         ratio = (preRed == 0) ? 0 : (actRed / preRed);
353 
354         // update the step bound
355         if (ratio <= 0.25) {
356           double tmp =
357             (actRed < 0) ? (0.5 * dirDer / (dirDer + 0.5 * actRed)) : 0.5;
358           if ((0.1 * cost >= previousCost) || (tmp < 0.1)) {
359             tmp = 0.1;
360           }
361           delta = tmp * Math.min(delta, 10.0 * lmNorm);
362           lmPar /= tmp;
363         } else if ((lmPar == 0) || (ratio >= 0.75)) {
364           delta = 2 * lmNorm;
365           lmPar *= 0.5;
366         }
367 
368         // test for successful iteration.
369         if (ratio >= 1.0e-4) {
370           // successful iteration, update the norm
371           firstIteration = false;
372           xNorm = 0;
373           for (int k = 0; k < cols; ++k) {
374             double xK = diag[k] * parameters[k].getEstimate();
375             xNorm    += xK * xK;
376           }
377           xNorm = Math.sqrt(xNorm);
378         } else {
379           // failed iteration, reset the previous values
380           cost = previousCost;
381           for (int j = 0; j < solvedCols; ++j) {
382             int pj = permutation[j];
383             parameters[pj].setEstimate(oldX[pj]);
384           }
385           tmpVec    = residuals;
386           residuals = oldRes;
387           oldRes    = tmpVec;
388         }
389    
390         // tests for convergence.
391         if (((Math.abs(actRed) <= costRelativeTolerance) &&
392              (preRed <= costRelativeTolerance) &&
393              (ratio <= 2.0)) ||
394              (delta <= parRelativeTolerance * xNorm)) {
395           return;
396         }
397 
398         // tests for termination and stringent tolerances
399         // (2.2204e-16 is the machine epsilon for IEEE754)
400         if ((Math.abs(actRed) <= 2.2204e-16) && (preRed <= 2.2204e-16) && (ratio <= 2.0)) {
401           throw new EstimationException("cost relative tolerance is too small ({0})," +
402                                         " no further reduction in the" +
403                                         " sum of squares is possible",
404                                         new Object[] { new Double(costRelativeTolerance) });
405         } else if (delta <= 2.2204e-16 * xNorm) {
406           throw new EstimationException("parameters relative tolerance is too small" +
407                                         " ({0}), no further improvement in" +
408                                         " the approximate solution is possible",
409                                         new Object[] { new Double(parRelativeTolerance) });
410         } else if (maxCosine <= 2.2204e-16)  {
411           throw new EstimationException("orthogonality tolerance is too small ({0})," +
412                                         " solution is orthogonal to the jacobian",
413                                         new Object[] { new Double(orthoTolerance) });
414         }
415 
416       }
417 
418     }
419 
420   }
421 
422   /** 
423    * Determine the Levenberg-Marquardt parameter.
424    * <p>This implementation is a translation in Java of the MINPACK
425    * <a href="http://www.netlib.org/minpack/lmpar.f">lmpar</a>
426    * routine.</p>
427    * <p>This method sets the lmPar and lmDir attributes.</p>
428    * <p>The authors of the original fortran function are:</p>
429    * <ul>
430    *   <li>Argonne National Laboratory. MINPACK project. March 1980</li>
431    *   <li>Burton  S. Garbow</li>
432    *   <li>Kenneth E. Hillstrom</li>
433    *   <li>Jorge   J. More</li>
434    * </ul>
435    * <p>Luc Maisonobe did the Java translation.</p>
436    * 
437    * @param qy array containing qTy
438    * @param delta upper bound on the euclidean norm of diagR * lmDir
439    * @param diag diagonal matrix
440    * @param work1 work array
441    * @param work2 work array
442    * @param work3 work array
443    */
444   private void determineLMParameter(double[] qy, double delta, double[] diag,
445                                     double[] work1, double[] work2, double[] work3) {
446 
447     // compute and store in x the gauss-newton direction, if the
448     // jacobian is rank-deficient, obtain a least squares solution
449     for (int j = 0; j < rank; ++j) {
450       lmDir[permutation[j]] = qy[j];
451     }
452     for (int j = rank; j < cols; ++j) {
453       lmDir[permutation[j]] = 0;
454     }
455     for (int k = rank - 1; k >= 0; --k) {
456       int pk = permutation[k];
457       double ypk = lmDir[pk] / diagR[pk];
458       for (int i = 0, index = pk; i < k; ++i, index += cols) {
459         lmDir[permutation[i]] -= ypk * jacobian[index];
460       }
461       lmDir[pk] = ypk;
462     }
463 
464     // evaluate the function at the origin, and test
465     // for acceptance of the Gauss-Newton direction
466     double dxNorm = 0;
467     for (int j = 0; j < solvedCols; ++j) {
468       int pj = permutation[j];
469       double s = diag[pj] * lmDir[pj];
470       work1[pj] = s;
471       dxNorm += s * s;
472     }
473     dxNorm = Math.sqrt(dxNorm);
474     double fp = dxNorm - delta;
475     if (fp <= 0.1 * delta) {
476       lmPar = 0;
477       return;
478     }
479 
480     // if the jacobian is not rank deficient, the Newton step provides
481     // a lower bound, parl, for the zero of the function,
482     // otherwise set this bound to zero
483     double sum2, parl = 0;
484     if (rank == solvedCols) {
485       for (int j = 0; j < solvedCols; ++j) {
486         int pj = permutation[j];
487         work1[pj] *= diag[pj] / dxNorm; 
488       }
489       sum2 = 0;
490       for (int j = 0; j < solvedCols; ++j) {
491         int pj = permutation[j];
492         double sum = 0;
493         for (int i = 0, index = pj; i < j; ++i, index += cols) {
494           sum += jacobian[index] * work1[permutation[i]];
495         }
496         double s = (work1[pj] - sum) / diagR[pj];
497         work1[pj] = s;
498         sum2 += s * s;
499       }
500       parl = fp / (delta * sum2);
501     }
502 
503     // calculate an upper bound, paru, for the zero of the function
504     sum2 = 0;
505     for (int j = 0; j < solvedCols; ++j) {
506       int pj = permutation[j];
507       double sum = 0;
508       for (int i = 0, index = pj; i <= j; ++i, index += cols) {
509         sum += jacobian[index] * qy[i];
510       }
511       sum /= diag[pj];
512       sum2 += sum * sum;
513     }
514     double gNorm = Math.sqrt(sum2);
515     double paru = gNorm / delta;
516     if (paru == 0) {
517       // 2.2251e-308 is the smallest positive real for IEE754
518       paru = 2.2251e-308 / Math.min(delta, 0.1);
519     }
520 
521     // if the input par lies outside of the interval (parl,paru),
522     // set par to the closer endpoint
523     lmPar = Math.min(paru, Math.max(lmPar, parl));
524     if (lmPar == 0) {
525       lmPar = gNorm / dxNorm;
526     }
527 
528     for (int countdown = 10; countdown >= 0; --countdown) {
529 
530       // evaluate the function at the current value of lmPar
531       if (lmPar == 0) {
532         lmPar = Math.max(2.2251e-308, 0.001 * paru);
533       }
534       double sPar = Math.sqrt(lmPar);
535       for (int j = 0; j < solvedCols; ++j) {
536         int pj = permutation[j];
537         work1[pj] = sPar * diag[pj];
538       }
539       determineLMDirection(qy, work1, work2, work3);
540 
541       dxNorm = 0;
542       for (int j = 0; j < solvedCols; ++j) {
543         int pj = permutation[j];
544         double s = diag[pj] * lmDir[pj];
545         work3[pj] = s;
546         dxNorm += s * s;
547       }
548       dxNorm = Math.sqrt(dxNorm);
549       double previousFP = fp;
550       fp = dxNorm - delta;
551 
552       // if the function is small enough, accept the current value
553       // of lmPar, also test for the exceptional cases where parl is zero
554       if ((Math.abs(fp) <= 0.1 * delta) ||
555           ((parl == 0) && (fp <= previousFP) && (previousFP < 0))) {
556         return;
557       }
558  
559       // compute the Newton correction
560       for (int j = 0; j < solvedCols; ++j) {
561        int pj = permutation[j];
562         work1[pj] = work3[pj] * diag[pj] / dxNorm; 
563       }
564       for (int j = 0; j < solvedCols; ++j) {
565         int pj = permutation[j];
566         work1[pj] /= work2[j];
567         double tmp = work1[pj];
568         for (int i = j + 1; i < solvedCols; ++i) {
569           work1[permutation[i]] -= jacobian[i * cols + pj] * tmp;
570         }
571       }
572       sum2 = 0;
573       for (int j = 0; j < solvedCols; ++j) {
574         double s = work1[permutation[j]];
575         sum2 += s * s;
576       }
577       double correction = fp / (delta * sum2);
578 
579       // depending on the sign of the function, update parl or paru.
580       if (fp > 0) {
581         parl = Math.max(parl, lmPar);
582       } else if (fp < 0) {
583         paru = Math.min(paru, lmPar);
584       }
585 
586       // compute an improved estimate for lmPar
587       lmPar = Math.max(parl, lmPar + correction);
588 
589     }
590   }
591 
592   /** 
593    * Solve a*x = b and d*x = 0 in the least squares sense.
594    * <p>This implementation is a translation in Java of the MINPACK
595    * <a href="http://www.netlib.org/minpack/qrsolv.f">qrsolv</a>
596    * routine.</p>
597    * <p>This method sets the lmDir and lmDiag attributes.</p>
598    * <p>The authors of the original fortran function are:</p>
599    * <ul>
600    *   <li>Argonne National Laboratory. MINPACK project. March 1980</li>
601    *   <li>Burton  S. Garbow</li>
602    *   <li>Kenneth E. Hillstrom</li>
603    *   <li>Jorge   J. More</li>
604    * </ul>
605    * <p>Luc Maisonobe did the Java translation.</p>
606    * 
607    * @param qy array containing qTy
608    * @param diag diagonal matrix
609    * @param lmDiag diagonal elements associated with lmDir
610    * @param work work array
611    */
612   private void determineLMDirection(double[] qy, double[] diag,
613                                     double[] lmDiag, double[] work) {
614 
615     // copy R and Qty to preserve input and initialize s
616     //  in particular, save the diagonal elements of R in lmDir
617     for (int j = 0; j < solvedCols; ++j) {
618       int pj = permutation[j];
619       for (int i = j + 1; i < solvedCols; ++i) {
620         jacobian[i * cols + pj] = jacobian[j * cols + permutation[i]];
621       }
622       lmDir[j] = diagR[pj];
623       work[j]  = qy[j];
624     }
625 
626     // eliminate the diagonal matrix d using a Givens rotation
627     for (int j = 0; j < solvedCols; ++j) {
628 
629       // prepare the row of d to be eliminated, locating the
630       // diagonal element using p from the Q.R. factorization
631       int pj = permutation[j];
632       double dpj = diag[pj];
633       if (dpj != 0) {
634         Arrays.fill(lmDiag, j + 1, lmDiag.length, 0);
635       }
636       lmDiag[j] = dpj;
637 
638       //  the transformations to eliminate the row of d
639       // modify only a single element of Qty
640       // beyond the first n, which is initially zero.
641       double qtbpj = 0;
642       for (int k = j; k < solvedCols; ++k) {
643         int pk = permutation[k];
644 
645         // determine a Givens rotation which eliminates the
646         // appropriate element in the current row of d
647         if (lmDiag[k] != 0) {
648 
649           double sin, cos;
650           double rkk = jacobian[k * cols + pk];
651           if (Math.abs(rkk) < Math.abs(lmDiag[k])) {
652             double cotan = rkk / lmDiag[k];
653             sin   = 1.0 / Math.sqrt(1.0 + cotan * cotan);
654             cos   = sin * cotan;
655           } else {
656             double tan = lmDiag[k] / rkk;
657             cos = 1.0 / Math.sqrt(1.0 + tan * tan);
658             sin = cos * tan;
659           }
660 
661           // compute the modified diagonal element of R and
662           // the modified element of (Qty,0)
663           jacobian[k * cols + pk] = cos * rkk + sin * lmDiag[k];
664           double temp = cos * work[k] + sin * qtbpj;
665           qtbpj = -sin * work[k] + cos * qtbpj;
666           work[k] = temp;
667 
668           // accumulate the tranformation in the row of s
669           for (int i = k + 1; i < solvedCols; ++i) {
670             double rik = jacobian[i * cols + pk];
671             temp = cos * rik + sin * lmDiag[i];
672             lmDiag[i] = -sin * rik + cos * lmDiag[i];
673             jacobian[i * cols + pk] = temp;
674           }
675 
676         }
677       }
678 
679       // store the diagonal element of s and restore
680       // the corresponding diagonal element of R
681       int index = j * cols + permutation[j];
682       lmDiag[j]       = jacobian[index];
683       jacobian[index] = lmDir[j];
684 
685     }
686 
687     // solve the triangular system for z, if the system is
688     // singular, then obtain a least squares solution
689     int nSing = solvedCols;
690     for (int j = 0; j < solvedCols; ++j) {
691       if ((lmDiag[j] == 0) && (nSing == solvedCols)) {
692         nSing = j;
693       }
694       if (nSing < solvedCols) {
695         work[j] = 0;
696       }
697     }
698     if (nSing > 0) {
699       for (int j = nSing - 1; j >= 0; --j) {
700         int pj = permutation[j];
701         double sum = 0;
702         for (int i = j + 1; i < nSing; ++i) {
703           sum += jacobian[i * cols + pj] * work[i];
704         }
705         work[j] = (work[j] - sum) / lmDiag[j];
706       }
707     }
708 
709     // permute the components of z back to components of lmDir
710     for (int j = 0; j < lmDir.length; ++j) {
711       lmDir[permutation[j]] = work[j];
712     }
713 
714   }
715 
716   /** 
717    * Decompose a matrix A as A.P = Q.R using Householder transforms.
718    * <p>As suggested in the P. Lascaux and R. Theodor book
719    * <i>Analyse num&eacute;rique matricielle appliqu&eacute;e &agrave;
720    * l'art de l'ing&eacute;nieur</i> (Masson, 1986), instead of representing
721    * the Householder transforms with u<sub>k</sub> unit vectors such that:
722    * <pre>
723    * H<sub>k</sub> = I - 2u<sub>k</sub>.u<sub>k</sub><sup>t</sup>
724    * </pre>
725    * we use <sub>k</sub> non-unit vectors such that:
726    * <pre>
727    * H<sub>k</sub> = I - beta<sub>k</sub>v<sub>k</sub>.v<sub>k</sub><sup>t</sup>
728    * </pre>
729    * where v<sub>k</sub> = a<sub>k</sub> - alpha<sub>k</sub> e<sub>k</sub>.
730    * The beta<sub>k</sub> coefficients are provided upon exit as recomputing
731    * them from the v<sub>k</sub> vectors would be costly.</p>
732    * <p>This decomposition handles rank deficient cases since the tranformations
733    * are performed in non-increasing columns norms order thanks to columns
734    * pivoting. The diagonal elements of the R matrix are therefore also in
735    * non-increasing absolute values order.</p>
736    */
737   private void qrDecomposition() {
738 
739     // initializations
740     for (int k = 0; k < cols; ++k) {
741       permutation[k] = k;
742       double norm2 = 0;
743       for (int index = k; index < jacobian.length; index += cols) {
744         double akk = jacobian[index];
745         norm2 += akk * akk;
746       }
747       jacNorm[k] = Math.sqrt(norm2);
748     }
749 
750     // transform the matrix column after column
751     for (int k = 0; k < cols; ++k) {
752 
753       // select the column with the greatest norm on active components
754       int nextColumn = -1;
755       double ak2 = Double.NEGATIVE_INFINITY;
756       for (int i = k; i < cols; ++i) {
757         double norm2 = 0;
758         int iDiag = k * cols + permutation[i];
759         for (int index = iDiag; index < jacobian.length; index += cols) {
760           double aki = jacobian[index];
761           norm2 += aki * aki;
762         }
763         if (norm2 > ak2) {
764           nextColumn = i;
765           ak2        = norm2;
766         }
767       }
768       if (ak2 == 0) {
769         rank = k;
770         return;
771       }
772       int pk                  = permutation[nextColumn];
773       permutation[nextColumn] = permutation[k];
774       permutation[k]          = pk;
775 
776       // choose alpha such that Hk.u = alpha ek
777       int    kDiag = k * cols + pk;
778       double akk   = jacobian[kDiag];
779       double alpha = (akk > 0) ? -Math.sqrt(ak2) : Math.sqrt(ak2);
780       double betak = 1.0 / (ak2 - akk * alpha);
781       beta[pk]     = betak;
782 
783       // transform the current column
784       diagR[pk]        = alpha;
785       jacobian[kDiag] -= alpha;
786 
787       // transform the remaining columns
788       for (int dk = cols - 1 - k; dk > 0; --dk) {
789         int dkp = permutation[k + dk] - pk;
790         double gamma = 0;
791         for (int index = kDiag; index < jacobian.length; index += cols) {
792           gamma += jacobian[index] * jacobian[index + dkp];
793         }
794         gamma *= betak;
795         for (int index = kDiag; index < jacobian.length; index += cols) {
796           jacobian[index + dkp] -= gamma * jacobian[index];
797         }
798       }
799 
800     }
801 
802     rank = solvedCols;
803 
804   }
805 
806   /** 
807    * Compute the product Qt.y for some Q.R. decomposition.
808    * 
809    * @param y vector to multiply (will be overwritten with the result)
810    */
811   private void qTy(double[] y) {
812     for (int k = 0; k < cols; ++k) {
813       int pk = permutation[k];
814       int kDiag = k * cols + pk;
815       double gamma = 0;
816       for (int i = k, index = kDiag; i < rows; ++i, index += cols) {
817         gamma += jacobian[index] * y[i];
818       }
819       gamma *= beta[pk];
820       for (int i = k, index = kDiag; i < rows; ++i, index += cols) {
821         y[i] -= gamma * jacobian[index];
822       }
823     }
824   }
825 
826   /** Number of solved variables. */
827   private int solvedCols;
828 
829   /** Diagonal elements of the R matrix in the Q.R. decomposition. */
830   private double[] diagR;
831 
832   /** Norms of the columns of the jacobian matrix. */
833   private double[] jacNorm;
834 
835   /** Coefficients of the Householder transforms vectors. */
836   private double[] beta;
837 
838   /** Columns permutation array. */
839   private int[] permutation;
840 
841   /** Rank of the jacobian matrix. */
842   private int rank;
843 
844   /** Levenberg-Marquardt parameter. */
845   private double lmPar;
846 
847   /** Parameters evolution direction associated with lmPar. */
848   private double[] lmDir;
849 
850   /** Positive input variable used in determining the initial step bound. */
851   private double initialStepBoundFactor;
852 
853   /** Desired relative error in the sum of squares. */
854   private double costRelativeTolerance;
855 
856   /**  Desired relative error in the approximate solution parameters. */
857   private double parRelativeTolerance;
858 
859   /** Desired max cosine on the orthogonality between the function vector
860    * and the columns of the jacobian. */
861   private double orthoTolerance;
862 
863   /** Serializable version identifier */
864   private static final long serialVersionUID = -5705952631533171019L;
865 
866 }