LeastSquaresFactory.java

  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.  *      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. /*
  18.  * This is not the original file distributed by the Apache Software Foundation
  19.  * It has been modified by the Hipparchus project
  20.  */
  21. package org.hipparchus.optim.nonlinear.vector.leastsquares;

  22. import org.hipparchus.analysis.MultivariateMatrixFunction;
  23. import org.hipparchus.analysis.MultivariateVectorFunction;
  24. import org.hipparchus.exception.MathIllegalStateException;
  25. import org.hipparchus.linear.Array2DRowRealMatrix;
  26. import org.hipparchus.linear.ArrayRealVector;
  27. import org.hipparchus.linear.DiagonalMatrix;
  28. import org.hipparchus.linear.EigenDecompositionSymmetric;
  29. import org.hipparchus.linear.RealMatrix;
  30. import org.hipparchus.linear.RealVector;
  31. import org.hipparchus.optim.AbstractOptimizationProblem;
  32. import org.hipparchus.optim.ConvergenceChecker;
  33. import org.hipparchus.optim.LocalizedOptimFormats;
  34. import org.hipparchus.optim.PointVectorValuePair;
  35. import org.hipparchus.optim.nonlinear.vector.leastsquares.LeastSquaresProblem.Evaluation;
  36. import org.hipparchus.util.FastMath;
  37. import org.hipparchus.util.Incrementor;
  38. import org.hipparchus.util.Pair;

  39. /**
  40.  * A Factory for creating {@link LeastSquaresProblem}s.
  41.  *
  42.  */
  43. public class LeastSquaresFactory {

  44.     /** Prevent instantiation. */
  45.     private LeastSquaresFactory() {}

  46.     /**
  47.      * Create a {@link org.hipparchus.optim.nonlinear.vector.leastsquares.LeastSquaresProblem}
  48.      * from the given elements. There will be no weights applied (unit weights).
  49.      *
  50.      * @param model          the model function. Produces the computed values.
  51.      * @param observed       the observed (target) values
  52.      * @param start          the initial guess.
  53.      * @param weight         the weight matrix
  54.      * @param checker        convergence checker
  55.      * @param maxEvaluations the maximum number of times to evaluate the model
  56.      * @param maxIterations  the maximum number to times to iterate in the algorithm
  57.      * @param lazyEvaluation Whether the call to {@link Evaluation#evaluate(RealVector)}
  58.      * will defer the evaluation until access to the value is requested.
  59.      * @param paramValidator Model parameters validator.
  60.      * @return the specified General Least Squares problem.
  61.      *
  62.      */
  63.     public static LeastSquaresProblem create(final MultivariateJacobianFunction model,
  64.                                              final RealVector observed,
  65.                                              final RealVector start,
  66.                                              final RealMatrix weight,
  67.                                              final ConvergenceChecker<Evaluation> checker,
  68.                                              final int maxEvaluations,
  69.                                              final int maxIterations,
  70.                                              final boolean lazyEvaluation,
  71.                                              final ParameterValidator paramValidator) {
  72.         final LeastSquaresProblem p = new LocalLeastSquaresProblem(model,
  73.                                                                    observed,
  74.                                                                    start,
  75.                                                                    checker,
  76.                                                                    maxEvaluations,
  77.                                                                    maxIterations,
  78.                                                                    lazyEvaluation,
  79.                                                                    paramValidator);
  80.         if (weight != null) {
  81.             return weightMatrix(p, weight);
  82.         } else {
  83.             return p;
  84.         }
  85.     }

  86.     /**
  87.      * Create a {@link org.hipparchus.optim.nonlinear.vector.leastsquares.LeastSquaresProblem}
  88.      * from the given elements. There will be no weights applied (unit weights).
  89.      *
  90.      * @param model          the model function. Produces the computed values.
  91.      * @param observed       the observed (target) values
  92.      * @param start          the initial guess.
  93.      * @param checker        convergence checker
  94.      * @param maxEvaluations the maximum number of times to evaluate the model
  95.      * @param maxIterations  the maximum number to times to iterate in the algorithm
  96.      * @return the specified General Least Squares problem.
  97.      */
  98.     public static LeastSquaresProblem create(final MultivariateJacobianFunction model,
  99.                                              final RealVector observed,
  100.                                              final RealVector start,
  101.                                              final ConvergenceChecker<Evaluation> checker,
  102.                                              final int maxEvaluations,
  103.                                              final int maxIterations) {
  104.         return create(model,
  105.                       observed,
  106.                       start,
  107.                       null,
  108.                       checker,
  109.                       maxEvaluations,
  110.                       maxIterations,
  111.                       false,
  112.                       null);
  113.     }

  114.     /**
  115.      * Create a {@link org.hipparchus.optim.nonlinear.vector.leastsquares.LeastSquaresProblem}
  116.      * from the given elements.
  117.      *
  118.      * @param model          the model function. Produces the computed values.
  119.      * @param observed       the observed (target) values
  120.      * @param start          the initial guess.
  121.      * @param weight         the weight matrix
  122.      * @param checker        convergence checker
  123.      * @param maxEvaluations the maximum number of times to evaluate the model
  124.      * @param maxIterations  the maximum number to times to iterate in the algorithm
  125.      * @return the specified General Least Squares problem.
  126.      */
  127.     public static LeastSquaresProblem create(final MultivariateJacobianFunction model,
  128.                                              final RealVector observed,
  129.                                              final RealVector start,
  130.                                              final RealMatrix weight,
  131.                                              final ConvergenceChecker<Evaluation> checker,
  132.                                              final int maxEvaluations,
  133.                                              final int maxIterations) {
  134.         return weightMatrix(create(model,
  135.                                    observed,
  136.                                    start,
  137.                                    checker,
  138.                                    maxEvaluations,
  139.                                    maxIterations),
  140.                             weight);
  141.     }

  142.     /**
  143.      * Create a {@link org.hipparchus.optim.nonlinear.vector.leastsquares.LeastSquaresProblem}
  144.      * from the given elements.
  145.      * <p>
  146.      * This factory method is provided for continuity with previous interfaces. Newer
  147.      * applications should use {@link #create(MultivariateJacobianFunction, RealVector,
  148.      * RealVector, ConvergenceChecker, int, int)}, or {@link #create(MultivariateJacobianFunction,
  149.      * RealVector, RealVector, RealMatrix, ConvergenceChecker, int, int)}.
  150.      *
  151.      * @param model          the model function. Produces the computed values.
  152.      * @param jacobian       the jacobian of the model with respect to the parameters
  153.      * @param observed       the observed (target) values
  154.      * @param start          the initial guess.
  155.      * @param weight         the weight matrix
  156.      * @param checker        convergence checker
  157.      * @param maxEvaluations the maximum number of times to evaluate the model
  158.      * @param maxIterations  the maximum number to times to iterate in the algorithm
  159.      * @return the specified General Least Squares problem.
  160.      */
  161.     public static LeastSquaresProblem create(final MultivariateVectorFunction model,
  162.                                              final MultivariateMatrixFunction jacobian,
  163.                                              final double[] observed,
  164.                                              final double[] start,
  165.                                              final RealMatrix weight,
  166.                                              final ConvergenceChecker<Evaluation> checker,
  167.                                              final int maxEvaluations,
  168.                                              final int maxIterations) {
  169.         return create(model(model, jacobian),
  170.                       new ArrayRealVector(observed, false),
  171.                       new ArrayRealVector(start, false),
  172.                       weight,
  173.                       checker,
  174.                       maxEvaluations,
  175.                       maxIterations);
  176.     }

  177.     /**
  178.      * Apply a dense weight matrix to the {@link LeastSquaresProblem}.
  179.      *
  180.      * @param problem the unweighted problem
  181.      * @param weights the matrix of weights
  182.      * @return a new {@link LeastSquaresProblem} with the weights applied. The original
  183.      *         {@code problem} is not modified.
  184.      */
  185.     public static LeastSquaresProblem weightMatrix(final LeastSquaresProblem problem,
  186.                                                    final RealMatrix weights) {
  187.         final RealMatrix weightSquareRoot = squareRoot(weights);
  188.         return new LeastSquaresAdapter(problem) {
  189.             /** {@inheritDoc} */
  190.             @Override
  191.             public Evaluation evaluate(final RealVector point) {
  192.                 return new DenseWeightedEvaluation(super.evaluate(point), weightSquareRoot);
  193.             }
  194.         };
  195.     }

  196.     /**
  197.      * Apply a diagonal weight matrix to the {@link LeastSquaresProblem}.
  198.      *
  199.      * @param problem the unweighted problem
  200.      * @param weights the diagonal of the weight matrix
  201.      * @return a new {@link LeastSquaresProblem} with the weights applied. The original
  202.      *         {@code problem} is not modified.
  203.      */
  204.     public static LeastSquaresProblem weightDiagonal(final LeastSquaresProblem problem,
  205.                                                      final RealVector weights) {
  206.         // TODO more efficient implementation
  207.         return weightMatrix(problem, new DiagonalMatrix(weights.toArray()));
  208.     }

  209.     /**
  210.      * Count the evaluations of a particular problem. The {@code counter} will be
  211.      * incremented every time {@link LeastSquaresProblem#evaluate(RealVector)} is called on
  212.      * the <em>returned</em> problem.
  213.      *
  214.      * @param problem the problem to track.
  215.      * @param counter the counter to increment.
  216.      * @return a least squares problem that tracks evaluations
  217.      */
  218.     public static LeastSquaresProblem countEvaluations(final LeastSquaresProblem problem,
  219.                                                        final Incrementor counter) {
  220.         return new LeastSquaresAdapter(problem) {

  221.             /** {@inheritDoc} */
  222.             @Override
  223.             public Evaluation evaluate(final RealVector point) {
  224.                 counter.increment();
  225.                 return super.evaluate(point);
  226.             }

  227.             // Delegate the rest.
  228.         };
  229.     }

  230.     /**
  231.      * View a convergence checker specified for a {@link PointVectorValuePair} as one
  232.      * specified for an {@link Evaluation}.
  233.      *
  234.      * @param checker the convergence checker to adapt.
  235.      * @return a convergence checker that delegates to {@code checker}.
  236.      */
  237.     public static ConvergenceChecker<Evaluation> evaluationChecker(final ConvergenceChecker<PointVectorValuePair> checker) {
  238.         return new ConvergenceChecker<Evaluation>() {
  239.             /** {@inheritDoc} */
  240.             @Override
  241.             public boolean converged(final int iteration,
  242.                                      final Evaluation previous,
  243.                                      final Evaluation current) {
  244.                 return checker.converged(
  245.                         iteration,
  246.                         new PointVectorValuePair(
  247.                                 previous.getPoint().toArray(),
  248.                                 previous.getResiduals().toArray(),
  249.                                 false),
  250.                         new PointVectorValuePair(
  251.                                 current.getPoint().toArray(),
  252.                                 current.getResiduals().toArray(),
  253.                                 false)
  254.                 );
  255.             }
  256.         };
  257.     }

  258.     /**
  259.      * Computes the square-root of the weight matrix.
  260.      *
  261.      * @param m Symmetric, positive-definite (weight) matrix.
  262.      * @return the square-root of the weight matrix.
  263.      */
  264.     private static RealMatrix squareRoot(final RealMatrix m) {
  265.         if (m instanceof DiagonalMatrix) {
  266.             final int dim = m.getRowDimension();
  267.             final RealMatrix sqrtM = new DiagonalMatrix(dim);
  268.             for (int i = 0; i < dim; i++) {
  269.                 sqrtM.setEntry(i, i, FastMath.sqrt(m.getEntry(i, i)));
  270.             }
  271.             return sqrtM;
  272.         } else {
  273.             final EigenDecompositionSymmetric dec = new EigenDecompositionSymmetric(m);
  274.             return dec.getSquareRoot();
  275.         }
  276.     }

  277.     /**
  278.      * Combine a {@link MultivariateVectorFunction} with a {@link
  279.      * MultivariateMatrixFunction} to produce a {@link MultivariateJacobianFunction}.
  280.      *
  281.      * @param value    the vector value function
  282.      * @param jacobian the Jacobian function
  283.      * @return a function that computes both at the same time
  284.      */
  285.     public static MultivariateJacobianFunction model(final MultivariateVectorFunction value,
  286.                                                      final MultivariateMatrixFunction jacobian) {
  287.         return new LocalValueAndJacobianFunction(value, jacobian);
  288.     }

  289.     /**
  290.      * Combine a {@link MultivariateVectorFunction} with a {@link
  291.      * MultivariateMatrixFunction} to produce a {@link MultivariateJacobianFunction}.
  292.      */
  293.     private static class LocalValueAndJacobianFunction
  294.         implements ValueAndJacobianFunction {
  295.         /** Model. */
  296.         private final MultivariateVectorFunction value;
  297.         /** Model's Jacobian. */
  298.         private final MultivariateMatrixFunction jacobian;

  299.         /**
  300.          * @param value Model function.
  301.          * @param jacobian Model's Jacobian function.
  302.          */
  303.         LocalValueAndJacobianFunction(final MultivariateVectorFunction value,
  304.                                       final MultivariateMatrixFunction jacobian) {
  305.             this.value = value;
  306.             this.jacobian = jacobian;
  307.         }

  308.         /** {@inheritDoc} */
  309.         @Override
  310.         public Pair<RealVector, RealMatrix> value(final RealVector point) {
  311.             //TODO get array from RealVector without copying?
  312.             final double[] p = point.toArray();

  313.             // Evaluate.
  314.             return new Pair<>(computeValue(p),
  315.                     computeJacobian(p));
  316.         }

  317.         /** {@inheritDoc} */
  318.         @Override
  319.         public RealVector computeValue(final double[] params) {
  320.             return new ArrayRealVector(value.value(params), false);
  321.         }

  322.         /** {@inheritDoc} */
  323.         @Override
  324.         public RealMatrix computeJacobian(final double[] params) {
  325.             return new Array2DRowRealMatrix(jacobian.value(params), false);
  326.         }
  327.     }


  328.     /**
  329.      * A private, "field" immutable (not "real" immutable) implementation of {@link
  330.      * LeastSquaresProblem}.
  331.      */
  332.     private static class LocalLeastSquaresProblem
  333.             extends AbstractOptimizationProblem<Evaluation>
  334.             implements LeastSquaresProblem {

  335.         /** Target values for the model function at optimum. */
  336.         private final RealVector target;
  337.         /** Model function. */
  338.         private final MultivariateJacobianFunction model;
  339.         /** Initial guess. */
  340.         private final RealVector start;
  341.         /** Whether to use lazy evaluation. */
  342.         private final boolean lazyEvaluation;
  343.         /** Model parameters validator. */
  344.         private final ParameterValidator paramValidator;

  345.         /**
  346.          * Create a {@link LeastSquaresProblem} from the given data.
  347.          *
  348.          * @param model          the model function
  349.          * @param target         the observed data
  350.          * @param start          the initial guess
  351.          * @param checker        the convergence checker
  352.          * @param maxEvaluations the allowed evaluations
  353.          * @param maxIterations  the allowed iterations
  354.          * @param lazyEvaluation Whether the call to {@link Evaluation#evaluate(RealVector)}
  355.          * will defer the evaluation until access to the value is requested.
  356.          * @param paramValidator Model parameters validator.
  357.          */
  358.         LocalLeastSquaresProblem(final MultivariateJacobianFunction model,
  359.                                  final RealVector target,
  360.                                  final RealVector start,
  361.                                  final ConvergenceChecker<Evaluation> checker,
  362.                                  final int maxEvaluations,
  363.                                  final int maxIterations,
  364.                                  final boolean lazyEvaluation,
  365.                                  final ParameterValidator paramValidator) {
  366.             super(maxEvaluations, maxIterations, checker);
  367.             this.target = target;
  368.             this.model = model;
  369.             this.start = start;
  370.             this.lazyEvaluation = lazyEvaluation;
  371.             this.paramValidator = paramValidator;

  372.             if (lazyEvaluation &&
  373.                 !(model instanceof ValueAndJacobianFunction)) {
  374.                 // Lazy evaluation requires that value and Jacobian
  375.                 // can be computed separately.
  376.                 throw new MathIllegalStateException(LocalizedOptimFormats.INVALID_IMPLEMENTATION,
  377.                                                     model.getClass().getName());
  378.             }
  379.         }

  380.         /** {@inheritDoc} */
  381.         @Override
  382.         public int getObservationSize() {
  383.             return target.getDimension();
  384.         }

  385.         /** {@inheritDoc} */
  386.         @Override
  387.         public int getParameterSize() {
  388.             return start.getDimension();
  389.         }

  390.         /** {@inheritDoc} */
  391.         @Override
  392.         public RealVector getStart() {
  393.             return start == null ? null : start.copy();
  394.         }

  395.         /** {@inheritDoc} */
  396.         @Override
  397.         public Evaluation evaluate(final RealVector point) {
  398.             // Copy so optimizer can change point without changing our instance.
  399.             final RealVector p = paramValidator == null ?
  400.                 point.copy() :
  401.                 paramValidator.validate(point.copy());

  402.             if (lazyEvaluation) {
  403.                 return new LazyUnweightedEvaluation((ValueAndJacobianFunction) model,
  404.                                                     target,
  405.                                                     p);
  406.             } else {
  407.                 // Evaluate value and jacobian in one function call.
  408.                 final Pair<RealVector, RealMatrix> value = model.value(p);
  409.                 return new UnweightedEvaluation(value.getFirst(),
  410.                                                 value.getSecond(),
  411.                                                 target,
  412.                                                 p);
  413.             }
  414.         }

  415.         /**
  416.          * Container with the model evaluation at a particular point.
  417.          */
  418.         private static class UnweightedEvaluation extends AbstractEvaluation {
  419.             /** Point of evaluation. */
  420.             private final RealVector point;
  421.             /** Derivative at point. */
  422.             private final RealMatrix jacobian;
  423.             /** Computed residuals. */
  424.             private final RealVector residuals;

  425.             /**
  426.              * Create an {@link Evaluation} with no weights.
  427.              *
  428.              * @param values   the computed function values
  429.              * @param jacobian the computed function Jacobian
  430.              * @param target   the observed values
  431.              * @param point    the abscissa
  432.              */
  433.             private UnweightedEvaluation(final RealVector values,
  434.                                          final RealMatrix jacobian,
  435.                                          final RealVector target,
  436.                                          final RealVector point) {
  437.                 super(target.getDimension());
  438.                 this.jacobian = jacobian;
  439.                 this.point = point;
  440.                 this.residuals = target.subtract(values);
  441.             }

  442.             /** {@inheritDoc} */
  443.             @Override
  444.             public RealMatrix getJacobian() {
  445.                 return jacobian;
  446.             }

  447.             /** {@inheritDoc} */
  448.             @Override
  449.             public RealVector getPoint() {
  450.                 return point;
  451.             }

  452.             /** {@inheritDoc} */
  453.             @Override
  454.             public RealVector getResiduals() {
  455.                 return residuals;
  456.             }
  457.         }

  458.         /**
  459.          * Container with the model <em>lazy</em> evaluation at a particular point.
  460.          */
  461.         private static class LazyUnweightedEvaluation extends AbstractEvaluation {
  462.             /** Point of evaluation. */
  463.             private final RealVector point;
  464.             /** Model and Jacobian functions. */
  465.             private final ValueAndJacobianFunction model;
  466.             /** Target values for the model function at optimum. */
  467.             private final RealVector target;

  468.             /**
  469.              * Create an {@link Evaluation} with no weights.
  470.              *
  471.              * @param model  the model function
  472.              * @param target the observed values
  473.              * @param point  the abscissa
  474.              */
  475.             private LazyUnweightedEvaluation(final ValueAndJacobianFunction model,
  476.                                              final RealVector target,
  477.                                              final RealVector point) {
  478.                 super(target.getDimension());
  479.                 // Safe to cast as long as we control usage of this class.
  480.                 this.model = model;
  481.                 this.point = point;
  482.                 this.target = target;
  483.             }

  484.             /** {@inheritDoc} */
  485.             @Override
  486.             public RealMatrix getJacobian() {
  487.                 return model.computeJacobian(point.toArray());
  488.             }

  489.             /** {@inheritDoc} */
  490.             @Override
  491.             public RealVector getPoint() {
  492.                 return point;
  493.             }

  494.             /** {@inheritDoc} */
  495.             @Override
  496.             public RealVector getResiduals() {
  497.                 return target.subtract(model.computeValue(point.toArray()));
  498.             }
  499.         }
  500.     }
  501. }