SequentialGaussNewtonOptimizer.java

  1. /*
  2.  * Licensed to the Hipparchus project 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 Hipparchus project 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.  *      https://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.hipparchus.optim.nonlinear.vector.leastsquares;

  18. import org.hipparchus.exception.LocalizedCoreFormats;
  19. import org.hipparchus.exception.MathIllegalArgumentException;
  20. import org.hipparchus.exception.MathIllegalStateException;
  21. import org.hipparchus.exception.NullArgumentException;
  22. import org.hipparchus.linear.ArrayRealVector;
  23. import org.hipparchus.linear.CholeskyDecomposition;
  24. import org.hipparchus.linear.MatrixDecomposer;
  25. import org.hipparchus.linear.MatrixUtils;
  26. import org.hipparchus.linear.QRDecomposer;
  27. import org.hipparchus.linear.RealMatrix;
  28. import org.hipparchus.linear.RealVector;
  29. import org.hipparchus.optim.ConvergenceChecker;
  30. import org.hipparchus.optim.LocalizedOptimFormats;
  31. import org.hipparchus.optim.nonlinear.vector.leastsquares.LeastSquaresProblem.Evaluation;
  32. import org.hipparchus.util.Incrementor;
  33. import org.hipparchus.util.Pair;

  34. /**
  35.  * Sequential Gauss-Newton least-squares solver.
  36.  * <p>
  37.  * This class solve a least-square problem by solving the normal equations of
  38.  * the linearized problem at each iteration.
  39.  * </p>
  40.  *
  41.  */
  42. public class SequentialGaussNewtonOptimizer implements LeastSquaresOptimizer {

  43.     /**
  44.      * The singularity threshold for matrix decompositions. Determines when a
  45.      * {@link MathIllegalStateException} is thrown. The current value was the
  46.      * default value for {@link org.hipparchus.linear.LUDecomposition}.
  47.      */
  48.     private static final double SINGULARITY_THRESHOLD = 1e-11;

  49.     /** Decomposer. */
  50.     private final MatrixDecomposer decomposer;

  51.     /** Indicates if normal equations should be formed explicitly. */
  52.     private final boolean formNormalEquations;

  53.     /** Old evaluation previously computed. */
  54.     private final Evaluation oldEvaluation;

  55.     /** Old jacobian previously computed. */
  56.     private final RealMatrix oldLhs;

  57.     /** Old residuals previously computed. */
  58.     private final RealVector oldRhs;

  59.     /**
  60.      * Create a sequential Gauss Newton optimizer.
  61.      * <p>
  62.      * The default for the algorithm is to use QR decomposition, not
  63.      * form normal equations and have no previous evaluation
  64.      * </p>
  65.      *
  66.      */
  67.     public SequentialGaussNewtonOptimizer() {
  68.         this(new QRDecomposer(SINGULARITY_THRESHOLD), false, null);
  69.     }

  70.     /**
  71.      * Create a sequential Gauss Newton optimizer that uses the given matrix
  72.      * decomposition algorithm to solve the normal equations.
  73.      * <p>
  74.      * The {@code decomposer} is used to solve J<sup>T</sup>Jx=J<sup>T</sup>r.
  75.      * </p>
  76.      *
  77.      * @param decomposer the decomposition algorithm to use.
  78.      * @param formNormalEquations whether the normal equations should be explicitly
  79.      *                            formed. If {@code true} then {@code decomposer} is used
  80.      *                            to solve J<sup>T</sup>Jx=J<sup>T</sup>r, otherwise
  81.      *                            {@code decomposer} is used to solve Jx=r. If {@code
  82.      *                            decomposer} can only solve square systems then this
  83.      *                            parameter should be {@code true}.
  84.      * @param evaluation old evaluation previously computed, null if there are no previous evaluations.
  85.      */
  86.     public SequentialGaussNewtonOptimizer(final MatrixDecomposer decomposer,
  87.                                           final boolean formNormalEquations,
  88.                                           final Evaluation evaluation) {
  89.         this.decomposer          = decomposer;
  90.         this.formNormalEquations = formNormalEquations;
  91.         this.oldEvaluation       = evaluation;
  92.         if (evaluation == null) {
  93.             this.oldLhs = null;
  94.             this.oldRhs = null;
  95.         } else {
  96.             if (formNormalEquations) {
  97.                 final Pair<RealMatrix, RealVector> normalEquation =
  98.                                 computeNormalMatrix(evaluation.getJacobian(), evaluation.getResiduals());
  99.                 // solve the linearized least squares problem
  100.                 this.oldLhs = normalEquation.getFirst();
  101.                 this.oldRhs = normalEquation.getSecond();
  102.             } else {
  103.                 this.oldLhs = evaluation.getJacobian();
  104.                 this.oldRhs = evaluation.getResiduals();
  105.             }
  106.         }
  107.     }

  108.     /**
  109.      * Get the matrix decomposition algorithm.
  110.      *
  111.      * @return the decomposition algorithm.
  112.      */
  113.     public MatrixDecomposer getDecomposer() {
  114.         return decomposer;
  115.     }

  116.     /**
  117.      * Configure the matrix decomposition algorithm.
  118.      *
  119.      * @param newDecomposer the decomposition algorithm to use.
  120.      * @return a new instance.
  121.      */
  122.     public SequentialGaussNewtonOptimizer withDecomposer(final MatrixDecomposer newDecomposer) {
  123.         return new SequentialGaussNewtonOptimizer(newDecomposer,
  124.                                                   this.isFormNormalEquations(),
  125.                                                   this.getOldEvaluation());
  126.     }

  127.     /**
  128.      * Get if the normal equations are explicitly formed.
  129.      *
  130.      * @return if the normal equations should be explicitly formed. If {@code true} then
  131.      * {@code decomposer} is used to solve J<sup>T</sup>Jx=J<sup>T</sup>r, otherwise
  132.      * {@code decomposer} is used to solve Jx=r.
  133.      */
  134.     public boolean isFormNormalEquations() {
  135.         return formNormalEquations;
  136.     }

  137.     /**
  138.      * Configure if the normal equations should be explicitly formed.
  139.      *
  140.      * @param newFormNormalEquations whether the normal equations should be explicitly
  141.      *                               formed. If {@code true} then {@code decomposer} is used
  142.      *                               to solve J<sup>T</sup>Jx=J<sup>T</sup>r, otherwise
  143.      *                               {@code decomposer} is used to solve Jx=r. If {@code
  144.      *                               decomposer} can only solve square systems then this
  145.      *                               parameter should be {@code true}.
  146.      * @return a new instance.
  147.      */
  148.     public SequentialGaussNewtonOptimizer withFormNormalEquations(final boolean newFormNormalEquations) {
  149.         return new SequentialGaussNewtonOptimizer(this.getDecomposer(),
  150.                                                   newFormNormalEquations,
  151.                                                   this.getOldEvaluation());
  152.     }

  153.     /**
  154.      * Get the previous evaluation used by the optimizer.
  155.      *
  156.      * @return the previous evaluation.
  157.      */
  158.     public Evaluation getOldEvaluation() {
  159.         return oldEvaluation;
  160.     }

  161.     /**
  162.      * Configure the previous evaluation used by the optimizer.
  163.      * <p>
  164.      * This building method uses a complete evaluation to retrieve
  165.      * a priori data. Note that as {@link #withAPrioriData(RealVector, RealMatrix)}
  166.      * generates a fake evaluation and calls this method, either
  167.      * {@link #withAPrioriData(RealVector, RealMatrix)} or {@link #withEvaluation(LeastSquaresProblem.Evaluation)}
  168.      * should be called, but not both as the last one called will override the previous one.
  169.      * </p>
  170.      * @param previousEvaluation the previous evaluation used by the optimizer.
  171.      * @return a new instance.
  172.      */
  173.     public SequentialGaussNewtonOptimizer withEvaluation(final Evaluation previousEvaluation) {
  174.         return new SequentialGaussNewtonOptimizer(this.getDecomposer(),
  175.                                                   this.isFormNormalEquations(),
  176.                                                   previousEvaluation);
  177.     }

  178.     /**
  179.      * Configure from a priori state and covariance.
  180.      * <p>
  181.      * This building method generates a fake evaluation and calls
  182.      * {@link #withEvaluation(LeastSquaresProblem.Evaluation)}, so either
  183.      * {@link #withAPrioriData(RealVector, RealMatrix)} or {@link #withEvaluation(LeastSquaresProblem.Evaluation)}
  184.      * should be called, but not both as the last one called will override the previous one.
  185.      * <p>
  186.      * A Cholesky decomposition is used to compute the weighted jacobian from the
  187.      * a priori covariance. This method uses the default thresholds of the decomposition.
  188.      * </p>
  189.      * @param aPrioriState a priori state to use
  190.      * @param aPrioriCovariance a priori covariance to use
  191.      * @return a new instance.
  192.      * @see #withAPrioriData(RealVector, RealMatrix, double, double)
  193.      */
  194.     public SequentialGaussNewtonOptimizer withAPrioriData(final RealVector aPrioriState,
  195.                                                           final RealMatrix aPrioriCovariance) {
  196.         return withAPrioriData(aPrioriState, aPrioriCovariance,
  197.                                CholeskyDecomposition.DEFAULT_RELATIVE_SYMMETRY_THRESHOLD,
  198.                                CholeskyDecomposition.DEFAULT_ABSOLUTE_POSITIVITY_THRESHOLD);
  199.     }

  200.     /**
  201.      * Configure from a priori state and covariance.
  202.      * <p>
  203.      * This building method generates a fake evaluation and calls
  204.      * {@link #withEvaluation(LeastSquaresProblem.Evaluation)}, so either
  205.      * {@link #withAPrioriData(RealVector, RealMatrix)} or {@link #withEvaluation(LeastSquaresProblem.Evaluation)}
  206.      * should be called, but not both as the last one called will override the previous one.
  207.      * <p>
  208.      * A Cholesky decomposition is used to compute the weighted jacobian from the
  209.      * a priori covariance.
  210.      * </p>
  211.      * @param aPrioriState a priori state to use
  212.      * @param aPrioriCovariance a priori covariance to use
  213.      * @param relativeSymmetryThreshold Cholesky decomposition threshold above which off-diagonal
  214.      * elements are considered too different and matrix not symmetric
  215.      * @param absolutePositivityThreshold Cholesky decomposition threshold below which diagonal
  216.      * elements are considered null and matrix not positive definite
  217.      * @return a new instance.
  218.      * @since 2.3
  219.      */
  220.     public SequentialGaussNewtonOptimizer withAPrioriData(final RealVector aPrioriState,
  221.                                                           final RealMatrix aPrioriCovariance,
  222.                                                           final double relativeSymmetryThreshold,
  223.                                                           final double absolutePositivityThreshold) {

  224.         // we consider the a priori state and covariance come from a
  225.         // previous estimation with exactly one observation of each state
  226.         // component, so partials are the identity matrix, weight is the
  227.         // square root of inverse of covariance, and residuals are zero

  228.         // create a fake weighted Jacobian
  229.         final RealMatrix jTj              = getDecomposer().decompose(aPrioriCovariance).getInverse();
  230.         final RealMatrix weightedJacobian = new CholeskyDecomposition(jTj,
  231.                                                                       relativeSymmetryThreshold,
  232.                                                                       absolutePositivityThreshold).getLT();

  233.         // create fake zero residuals
  234.         final RealVector residuals        = MatrixUtils.createRealVector(aPrioriState.getDimension());

  235.         // combine everything as an evaluation
  236.         final Evaluation fakeEvaluation   = new AbstractEvaluation(aPrioriState.getDimension()) {

  237.             /** {@inheritDoc} */
  238.             @Override
  239.             public RealVector getResiduals() {
  240.                 return residuals;
  241.             }

  242.             /** {@inheritDoc} */
  243.             @Override
  244.             public RealVector getPoint() {
  245.                 return aPrioriState;
  246.             }

  247.             /** {@inheritDoc} */
  248.             @Override
  249.             public RealMatrix getJacobian() {
  250.                 return weightedJacobian;
  251.             }
  252.         };

  253.         return withEvaluation(fakeEvaluation);

  254.     }

  255.     /** {@inheritDoc} */
  256.     @Override
  257.     public Optimum optimize(final LeastSquaresProblem lsp) {
  258.         // create local evaluation and iteration counts
  259.         final Incrementor evaluationCounter = lsp.getEvaluationCounter();
  260.         final Incrementor iterationCounter = lsp.getIterationCounter();
  261.         final ConvergenceChecker<Evaluation> checker =
  262.             lsp.getConvergenceChecker();

  263.         // Computation will be useless without a checker (see "for-loop").
  264.         if (checker == null) {
  265.             throw new NullArgumentException();
  266.         }

  267.         RealVector currentPoint = lsp.getStart();

  268.         if (oldEvaluation != null &&
  269.             currentPoint.getDimension() != oldEvaluation.getPoint().getDimension()) {
  270.             throw new MathIllegalStateException(LocalizedCoreFormats.DIMENSIONS_MISMATCH,
  271.                       currentPoint.getDimension(), oldEvaluation.getPoint().getDimension());
  272.         }

  273.         // iterate until convergence is reached
  274.         Evaluation current = null;
  275.         while (true) {
  276.             iterationCounter.increment();

  277.             // evaluate the objective function and its jacobian
  278.             final Evaluation previous = current;

  279.             // Value of the objective function at "currentPoint".
  280.             evaluationCounter.increment();
  281.             current = lsp.evaluate(currentPoint);
  282.             final RealVector currentResiduals = current.getResiduals();
  283.             final RealMatrix weightedJacobian = current.getJacobian();

  284.             currentPoint = current.getPoint();

  285.             // Check convergence.
  286.             if (previous != null &&
  287.                 checker.converged(iterationCounter.getCount(), previous,
  288.                                   current)) {
  289.                 // combine old and new evaluations
  290.                 final Evaluation combinedEvaluation = oldEvaluation == null ?
  291.                                                       current :
  292.                                                       new CombinedEvaluation(oldEvaluation, current);
  293.                 return Optimum.of(combinedEvaluation, evaluationCounter.getCount(),
  294.                                   iterationCounter.getCount());
  295.             }

  296.            // solve the linearized least squares problem
  297.             final RealMatrix lhs; // left hand side
  298.             final RealVector rhs; // right hand side
  299.             if (this.formNormalEquations) {
  300.                 final Pair<RealMatrix, RealVector> normalEquation =
  301.                                 computeNormalMatrix(weightedJacobian, currentResiduals);

  302.                 lhs = oldLhs == null ?
  303.                       normalEquation.getFirst() :
  304.                       normalEquation.getFirst().add(oldLhs); // left hand side
  305.                 rhs = oldRhs == null ?
  306.                       normalEquation.getSecond() :
  307.                       normalEquation.getSecond().add(oldRhs); // right hand side
  308.             } else {
  309.                 lhs = oldLhs == null ?
  310.                       weightedJacobian :
  311.                       combineJacobians(oldLhs, weightedJacobian);
  312.                 rhs = oldRhs == null ?
  313.                       currentResiduals :
  314.                       combineResiduals(oldRhs, currentResiduals);
  315.             }

  316.             final RealVector dX;
  317.             try {
  318.                 dX = this.decomposer.decompose(lhs).solve(rhs);
  319.             } catch (MathIllegalArgumentException e) {
  320.                 // change exception message
  321.                 throw new MathIllegalStateException(LocalizedOptimFormats.UNABLE_TO_SOLVE_SINGULAR_PROBLEM,
  322.                                                     e);
  323.             }
  324.             // update the estimated parameters
  325.             currentPoint = currentPoint.add(dX);

  326.         }
  327.     }

  328.     /** {@inheritDoc} */
  329.     @Override
  330.     public String toString() {
  331.         return "SequentialGaussNewtonOptimizer{" +
  332.                "decomposer=" + decomposer + '}';
  333.     }

  334.     /**
  335.      * Compute the normal matrix, J<sup>T</sup>J.
  336.      *
  337.      * @param jacobian  the m by n jacobian matrix, J. Input.
  338.      * @param residuals the m by 1 residual vector, r. Input.
  339.      * @return  the n by n normal matrix and the n by 1 J<sup>Tr</sup> vector.
  340.      */
  341.     private static Pair<RealMatrix, RealVector>
  342.         computeNormalMatrix(final RealMatrix jacobian,
  343.                             final RealVector residuals) {
  344.         // since the normal matrix is symmetric, we only need to compute half of
  345.         // it.
  346.         final int nR = jacobian.getRowDimension();
  347.         final int nC = jacobian.getColumnDimension();
  348.         // allocate space for return values
  349.         final RealMatrix normal = MatrixUtils.createRealMatrix(nC, nC);
  350.         final RealVector jTr = new ArrayRealVector(nC);
  351.         // for each measurement
  352.         for (int i = 0; i < nR; ++i) {
  353.             // compute JTr for measurement i
  354.             for (int j = 0; j < nC; j++) {
  355.                 jTr.setEntry(j,
  356.                              jTr.getEntry(j) +
  357.                                 residuals.getEntry(i) *
  358.                                                jacobian.getEntry(i, j));
  359.             }

  360.             // add the the contribution to the normal matrix for measurement i
  361.             for (int k = 0; k < nC; ++k) {
  362.                 // only compute the upper triangular part
  363.                 for (int l = k; l < nC; ++l) {
  364.                     normal
  365.                         .setEntry(k, l,
  366.                                   normal.getEntry(k,
  367.                                                   l) +
  368.                                         jacobian.getEntry(i, k) *
  369.                                                        jacobian.getEntry(i, l));
  370.                 }
  371.             }
  372.         }
  373.         // copy the upper triangular part to the lower triangular part.
  374.         for (int i = 0; i < nC; i++) {
  375.             for (int j = 0; j < i; j++) {
  376.                 normal.setEntry(i, j, normal.getEntry(j, i));
  377.             }
  378.         }
  379.         return new Pair<>(normal, jTr);
  380.     }

  381.     /** Combine Jacobian matrices
  382.      * @param oldJacobian old Jacobian matrix
  383.      * @param newJacobian new Jacobian matrix
  384.      * @return combined Jacobian matrix
  385.      */
  386.     private static RealMatrix combineJacobians(final RealMatrix oldJacobian,
  387.                                                final RealMatrix newJacobian) {
  388.         final int oldRowDimension    = oldJacobian.getRowDimension();
  389.         final int oldColumnDimension = oldJacobian.getColumnDimension();
  390.         final RealMatrix jacobian =
  391.                         MatrixUtils.createRealMatrix(oldRowDimension + newJacobian.getRowDimension(),
  392.                                                      oldColumnDimension);
  393.         jacobian.setSubMatrix(oldJacobian.getData(), 0,               0);
  394.         jacobian.setSubMatrix(newJacobian.getData(), oldRowDimension, 0);
  395.         return jacobian;
  396.     }

  397.     /** Combine residuals vectors
  398.      * @param oldResiduals old residuals vector
  399.      * @param newResiduals new residuals vector
  400.      * @return combined residuals vector
  401.      */
  402.     private static RealVector combineResiduals(final RealVector oldResiduals,
  403.                                                final RealVector newResiduals) {
  404.         return oldResiduals.append(newResiduals);
  405.     }

  406.     /**
  407.      * Container with an old and a new evaluation and combine both of them
  408.      */
  409.     private static class CombinedEvaluation extends AbstractEvaluation {

  410.         /** Point of evaluation. */
  411.         private final RealVector point;

  412.         /** Derivative at point. */
  413.         private final RealMatrix jacobian;

  414.         /** Computed residuals. */
  415.         private final RealVector residuals;

  416.         /**
  417.          * Create an {@link Evaluation} with no weights.
  418.          *
  419.          * @param oldEvaluation the old evaluation.
  420.          * @param newEvaluation the new evaluation
  421.          */
  422.         private CombinedEvaluation(final Evaluation oldEvaluation,
  423.                                    final Evaluation newEvaluation) {

  424.             super(oldEvaluation.getResiduals().getDimension() +
  425.                   newEvaluation.getResiduals().getDimension());

  426.             this.point    = newEvaluation.getPoint();
  427.             this.jacobian = combineJacobians(oldEvaluation.getJacobian(),
  428.                                              newEvaluation.getJacobian());
  429.             this.residuals = combineResiduals(oldEvaluation.getResiduals(),
  430.                                               newEvaluation.getResiduals());
  431.         }

  432.         /** {@inheritDoc} */
  433.         @Override
  434.         public RealMatrix getJacobian() {
  435.             return jacobian;
  436.         }

  437.         /** {@inheritDoc} */
  438.         @Override
  439.         public RealVector getPoint() {
  440.             return point;
  441.         }

  442.         /** {@inheritDoc} */
  443.         @Override
  444.         public RealVector getResiduals() {
  445.             return residuals;
  446.         }

  447.     }

  448. }