001/*
002 * Licensed to the Apache Software Foundation (ASF) under one or more
003 * contributor license agreements.  See the NOTICE file distributed with
004 * this work for additional information regarding copyright ownership.
005 * The ASF licenses this file to You under the Apache License, Version 2.0
006 * (the "License"); you may not use this file except in compliance with
007 * the License.  You may obtain a copy of the License at
008 *
009 *      http://www.apache.org/licenses/LICENSE-2.0
010 *
011 * Unless required by applicable law or agreed to in writing, software
012 * distributed under the License is distributed on an "AS IS" BASIS,
013 * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
014 * See the License for the specific language governing permissions and
015 * limitations under the License.
016 */
017package org.apache.commons.math3.fitting.leastsquares;
018
019import org.apache.commons.math3.exception.MathIllegalStateException;
020import org.apache.commons.math3.exception.util.LocalizedFormats;
021import org.apache.commons.math3.analysis.MultivariateMatrixFunction;
022import org.apache.commons.math3.analysis.MultivariateVectorFunction;
023import org.apache.commons.math3.fitting.leastsquares.LeastSquaresProblem.Evaluation;
024import org.apache.commons.math3.linear.Array2DRowRealMatrix;
025import org.apache.commons.math3.linear.ArrayRealVector;
026import org.apache.commons.math3.linear.DiagonalMatrix;
027import org.apache.commons.math3.linear.EigenDecomposition;
028import org.apache.commons.math3.linear.RealMatrix;
029import org.apache.commons.math3.linear.RealVector;
030import org.apache.commons.math3.optim.AbstractOptimizationProblem;
031import org.apache.commons.math3.optim.ConvergenceChecker;
032import org.apache.commons.math3.optim.PointVectorValuePair;
033import org.apache.commons.math3.util.FastMath;
034import org.apache.commons.math3.util.Incrementor;
035import org.apache.commons.math3.util.Pair;
036
037/**
038 * A Factory for creating {@link LeastSquaresProblem}s.
039 *
040 * @since 3.3
041 */
042public class LeastSquaresFactory {
043
044    /** Prevent instantiation. */
045    private LeastSquaresFactory() {}
046
047    /**
048     * Create a {@link org.apache.commons.math3.fitting.leastsquares.LeastSquaresProblem}
049     * from the given elements. There will be no weights applied (unit weights).
050     *
051     * @param model          the model function. Produces the computed values.
052     * @param observed       the observed (target) values
053     * @param start          the initial guess.
054     * @param weight         the weight matrix
055     * @param checker        convergence checker
056     * @param maxEvaluations the maximum number of times to evaluate the model
057     * @param maxIterations  the maximum number to times to iterate in the algorithm
058     * @param lazyEvaluation Whether the call to {@link Evaluation#evaluate(RealVector)}
059     * will defer the evaluation until access to the value is requested.
060     * @param paramValidator Model parameters validator.
061     * @return the specified General Least Squares problem.
062     *
063     * @since 3.4
064     */
065    public static LeastSquaresProblem create(final MultivariateJacobianFunction model,
066                                             final RealVector observed,
067                                             final RealVector start,
068                                             final RealMatrix weight,
069                                             final ConvergenceChecker<Evaluation> checker,
070                                             final int maxEvaluations,
071                                             final int maxIterations,
072                                             final boolean lazyEvaluation,
073                                             final ParameterValidator paramValidator) {
074        final LeastSquaresProblem p = new LocalLeastSquaresProblem(model,
075                                                                   observed,
076                                                                   start,
077                                                                   checker,
078                                                                   maxEvaluations,
079                                                                   maxIterations,
080                                                                   lazyEvaluation,
081                                                                   paramValidator);
082        if (weight != null) {
083            return weightMatrix(p, weight);
084        } else {
085            return p;
086        }
087    }
088
089    /**
090     * Create a {@link org.apache.commons.math3.fitting.leastsquares.LeastSquaresProblem}
091     * from the given elements. There will be no weights applied (unit weights).
092     *
093     * @param model          the model function. Produces the computed values.
094     * @param observed       the observed (target) values
095     * @param start          the initial guess.
096     * @param checker        convergence checker
097     * @param maxEvaluations the maximum number of times to evaluate the model
098     * @param maxIterations  the maximum number to times to iterate in the algorithm
099     * @return the specified General Least Squares problem.
100     */
101    public static LeastSquaresProblem create(final MultivariateJacobianFunction model,
102                                             final RealVector observed,
103                                             final RealVector start,
104                                             final ConvergenceChecker<Evaluation> checker,
105                                             final int maxEvaluations,
106                                             final int maxIterations) {
107        return create(model,
108                      observed,
109                      start,
110                      null,
111                      checker,
112                      maxEvaluations,
113                      maxIterations,
114                      false,
115                      null);
116    }
117
118    /**
119     * Create a {@link org.apache.commons.math3.fitting.leastsquares.LeastSquaresProblem}
120     * from the given elements.
121     *
122     * @param model          the model function. Produces the computed values.
123     * @param observed       the observed (target) values
124     * @param start          the initial guess.
125     * @param weight         the weight matrix
126     * @param checker        convergence checker
127     * @param maxEvaluations the maximum number of times to evaluate the model
128     * @param maxIterations  the maximum number to times to iterate in the algorithm
129     * @return the specified General Least Squares problem.
130     */
131    public static LeastSquaresProblem create(final MultivariateJacobianFunction model,
132                                             final RealVector observed,
133                                             final RealVector start,
134                                             final RealMatrix weight,
135                                             final ConvergenceChecker<Evaluation> checker,
136                                             final int maxEvaluations,
137                                             final int maxIterations) {
138        return weightMatrix(create(model,
139                                   observed,
140                                   start,
141                                   checker,
142                                   maxEvaluations,
143                                   maxIterations),
144                            weight);
145    }
146
147    /**
148     * Create a {@link org.apache.commons.math3.fitting.leastsquares.LeastSquaresProblem}
149     * from the given elements.
150     * <p>
151     * This factory method is provided for continuity with previous interfaces. Newer
152     * applications should use {@link #create(MultivariateJacobianFunction, RealVector,
153     * RealVector, ConvergenceChecker, int, int)}, or {@link #create(MultivariateJacobianFunction,
154     * RealVector, RealVector, RealMatrix, ConvergenceChecker, int, int)}.
155     *
156     * @param model          the model function. Produces the computed values.
157     * @param jacobian       the jacobian of the model with respect to the parameters
158     * @param observed       the observed (target) values
159     * @param start          the initial guess.
160     * @param weight         the weight matrix
161     * @param checker        convergence checker
162     * @param maxEvaluations the maximum number of times to evaluate the model
163     * @param maxIterations  the maximum number to times to iterate in the algorithm
164     * @return the specified General Least Squares problem.
165     */
166    public static LeastSquaresProblem create(final MultivariateVectorFunction model,
167                                             final MultivariateMatrixFunction jacobian,
168                                             final double[] observed,
169                                             final double[] start,
170                                             final RealMatrix weight,
171                                             final ConvergenceChecker<Evaluation> checker,
172                                             final int maxEvaluations,
173                                             final int maxIterations) {
174        return create(model(model, jacobian),
175                      new ArrayRealVector(observed, false),
176                      new ArrayRealVector(start, false),
177                      weight,
178                      checker,
179                      maxEvaluations,
180                      maxIterations);
181    }
182
183    /**
184     * Apply a dense weight matrix to the {@link LeastSquaresProblem}.
185     *
186     * @param problem the unweighted problem
187     * @param weights the matrix of weights
188     * @return a new {@link LeastSquaresProblem} with the weights applied. The original
189     *         {@code problem} is not modified.
190     */
191    public static LeastSquaresProblem weightMatrix(final LeastSquaresProblem problem,
192                                                   final RealMatrix weights) {
193        final RealMatrix weightSquareRoot = squareRoot(weights);
194        return new LeastSquaresAdapter(problem) {
195            /** {@inheritDoc} */
196            @Override
197            public Evaluation evaluate(final RealVector point) {
198                return new DenseWeightedEvaluation(super.evaluate(point), weightSquareRoot);
199            }
200        };
201    }
202
203    /**
204     * Apply a diagonal weight matrix to the {@link LeastSquaresProblem}.
205     *
206     * @param problem the unweighted problem
207     * @param weights the diagonal of the weight matrix
208     * @return a new {@link LeastSquaresProblem} with the weights applied. The original
209     *         {@code problem} is not modified.
210     */
211    public static LeastSquaresProblem weightDiagonal(final LeastSquaresProblem problem,
212                                                     final RealVector weights) {
213        // TODO more efficient implementation
214        return weightMatrix(problem, new DiagonalMatrix(weights.toArray()));
215    }
216
217    /**
218     * Count the evaluations of a particular problem. The {@code counter} will be
219     * incremented every time {@link LeastSquaresProblem#evaluate(RealVector)} is called on
220     * the <em>returned</em> problem.
221     *
222     * @param problem the problem to track.
223     * @param counter the counter to increment.
224     * @return a least squares problem that tracks evaluations
225     */
226    public static LeastSquaresProblem countEvaluations(final LeastSquaresProblem problem,
227                                                       final Incrementor counter) {
228        return new LeastSquaresAdapter(problem) {
229
230            /** {@inheritDoc} */
231            @Override
232            public Evaluation evaluate(final RealVector point) {
233                counter.incrementCount();
234                return super.evaluate(point);
235            }
236
237            // Delegate the rest.
238        };
239    }
240
241    /**
242     * View a convergence checker specified for a {@link PointVectorValuePair} as one
243     * specified for an {@link Evaluation}.
244     *
245     * @param checker the convergence checker to adapt.
246     * @return a convergence checker that delegates to {@code checker}.
247     */
248    public static ConvergenceChecker<Evaluation> evaluationChecker(final ConvergenceChecker<PointVectorValuePair> checker) {
249        return new ConvergenceChecker<Evaluation>() {
250            /** {@inheritDoc} */
251            public boolean converged(final int iteration,
252                                     final Evaluation previous,
253                                     final Evaluation current) {
254                return checker.converged(
255                        iteration,
256                        new PointVectorValuePair(
257                                previous.getPoint().toArray(),
258                                previous.getResiduals().toArray(),
259                                false),
260                        new PointVectorValuePair(
261                                current.getPoint().toArray(),
262                                current.getResiduals().toArray(),
263                                false)
264                );
265            }
266        };
267    }
268
269    /**
270     * Computes the square-root of the weight matrix.
271     *
272     * @param m Symmetric, positive-definite (weight) matrix.
273     * @return the square-root of the weight matrix.
274     */
275    private static RealMatrix squareRoot(final RealMatrix m) {
276        if (m instanceof DiagonalMatrix) {
277            final int dim = m.getRowDimension();
278            final RealMatrix sqrtM = new DiagonalMatrix(dim);
279            for (int i = 0; i < dim; i++) {
280                sqrtM.setEntry(i, i, FastMath.sqrt(m.getEntry(i, i)));
281            }
282            return sqrtM;
283        } else {
284            final EigenDecomposition dec = new EigenDecomposition(m);
285            return dec.getSquareRoot();
286        }
287    }
288
289    /**
290     * Combine a {@link MultivariateVectorFunction} with a {@link
291     * MultivariateMatrixFunction} to produce a {@link MultivariateJacobianFunction}.
292     *
293     * @param value    the vector value function
294     * @param jacobian the Jacobian function
295     * @return a function that computes both at the same time
296     */
297    public static MultivariateJacobianFunction model(final MultivariateVectorFunction value,
298                                                     final MultivariateMatrixFunction jacobian) {
299        return new LocalValueAndJacobianFunction(value, jacobian);
300    }
301
302    /**
303     * Combine a {@link MultivariateVectorFunction} with a {@link
304     * MultivariateMatrixFunction} to produce a {@link MultivariateJacobianFunction}.
305     *
306     * @param value    the vector value function
307     * @param jacobian the Jacobian function
308     * @return a function that computes both at the same time
309     */
310    private static class LocalValueAndJacobianFunction
311        implements ValueAndJacobianFunction {
312        /** Model. */
313        private final MultivariateVectorFunction value;
314        /** Model's Jacobian. */
315        private final MultivariateMatrixFunction jacobian;
316
317        /**
318         * @param value Model function.
319         * @param jacobian Model's Jacobian function.
320         */
321        LocalValueAndJacobianFunction(final MultivariateVectorFunction value,
322                                      final MultivariateMatrixFunction jacobian) {
323            this.value = value;
324            this.jacobian = jacobian;
325        }
326
327        /** {@inheritDoc} */
328        public Pair<RealVector, RealMatrix> value(final RealVector point) {
329            //TODO get array from RealVector without copying?
330            final double[] p = point.toArray();
331
332            // Evaluate.
333            return new Pair<RealVector, RealMatrix>(computeValue(p),
334                                                    computeJacobian(p));
335        }
336
337        /** {@inheritDoc} */
338        public RealVector computeValue(final double[] params) {
339            return new ArrayRealVector(value.value(params), false);
340        }
341
342        /** {@inheritDoc} */
343        public RealMatrix computeJacobian(final double[] params) {
344            return new Array2DRowRealMatrix(jacobian.value(params), false);
345        }
346    }
347
348
349    /**
350     * A private, "field" immutable (not "real" immutable) implementation of {@link
351     * LeastSquaresProblem}.
352     * @since 3.3
353     */
354    private static class LocalLeastSquaresProblem
355            extends AbstractOptimizationProblem<Evaluation>
356            implements LeastSquaresProblem {
357
358        /** Target values for the model function at optimum. */
359        private final RealVector target;
360        /** Model function. */
361        private final MultivariateJacobianFunction model;
362        /** Initial guess. */
363        private final RealVector start;
364        /** Whether to use lazy evaluation. */
365        private final boolean lazyEvaluation;
366        /** Model parameters validator. */
367        private final ParameterValidator paramValidator;
368
369        /**
370         * Create a {@link LeastSquaresProblem} from the given data.
371         *
372         * @param model          the model function
373         * @param target         the observed data
374         * @param start          the initial guess
375         * @param checker        the convergence checker
376         * @param maxEvaluations the allowed evaluations
377         * @param maxIterations  the allowed iterations
378         * @param lazyEvaluation Whether the call to {@link Evaluation#evaluate(RealVector)}
379         * will defer the evaluation until access to the value is requested.
380         * @param paramValidator Model parameters validator.
381         */
382        LocalLeastSquaresProblem(final MultivariateJacobianFunction model,
383                                 final RealVector target,
384                                 final RealVector start,
385                                 final ConvergenceChecker<Evaluation> checker,
386                                 final int maxEvaluations,
387                                 final int maxIterations,
388                                 final boolean lazyEvaluation,
389                                 final ParameterValidator paramValidator) {
390            super(maxEvaluations, maxIterations, checker);
391            this.target = target;
392            this.model = model;
393            this.start = start;
394            this.lazyEvaluation = lazyEvaluation;
395            this.paramValidator = paramValidator;
396
397            if (lazyEvaluation &&
398                !(model instanceof ValueAndJacobianFunction)) {
399                // Lazy evaluation requires that value and Jacobian
400                // can be computed separately.
401                throw new MathIllegalStateException(LocalizedFormats.INVALID_IMPLEMENTATION,
402                                                    model.getClass().getName());
403            }
404        }
405
406        /** {@inheritDoc} */
407        public int getObservationSize() {
408            return target.getDimension();
409        }
410
411        /** {@inheritDoc} */
412        public int getParameterSize() {
413            return start.getDimension();
414        }
415
416        /** {@inheritDoc} */
417        public RealVector getStart() {
418            return start == null ? null : start.copy();
419        }
420
421        /** {@inheritDoc} */
422        public Evaluation evaluate(final RealVector point) {
423            // Copy so optimizer can change point without changing our instance.
424            final RealVector p = paramValidator == null ?
425                point.copy() :
426                paramValidator.validate(point.copy());
427
428            if (lazyEvaluation) {
429                return new LazyUnweightedEvaluation((ValueAndJacobianFunction) model,
430                                                    target,
431                                                    p);
432            } else {
433                // Evaluate value and jacobian in one function call.
434                final Pair<RealVector, RealMatrix> value = model.value(p);
435                return new UnweightedEvaluation(value.getFirst(),
436                                                value.getSecond(),
437                                                target,
438                                                p);
439            }
440        }
441
442        /**
443         * Container with the model evaluation at a particular point.
444         */
445        private static class UnweightedEvaluation extends AbstractEvaluation {
446            /** Point of evaluation. */
447            private final RealVector point;
448            /** Derivative at point. */
449            private final RealMatrix jacobian;
450            /** Computed residuals. */
451            private final RealVector residuals;
452
453            /**
454             * Create an {@link Evaluation} with no weights.
455             *
456             * @param values   the computed function values
457             * @param jacobian the computed function Jacobian
458             * @param target   the observed values
459             * @param point    the abscissa
460             */
461            private UnweightedEvaluation(final RealVector values,
462                                         final RealMatrix jacobian,
463                                         final RealVector target,
464                                         final RealVector point) {
465                super(target.getDimension());
466                this.jacobian = jacobian;
467                this.point = point;
468                this.residuals = target.subtract(values);
469            }
470
471            /** {@inheritDoc} */
472            public RealMatrix getJacobian() {
473                return jacobian;
474            }
475
476            /** {@inheritDoc} */
477            public RealVector getPoint() {
478                return point;
479            }
480
481            /** {@inheritDoc} */
482            public RealVector getResiduals() {
483                return residuals;
484            }
485        }
486
487        /**
488         * Container with the model <em>lazy</em> evaluation at a particular point.
489         */
490        private static class LazyUnweightedEvaluation extends AbstractEvaluation {
491            /** Point of evaluation. */
492            private final RealVector point;
493            /** Model and Jacobian functions. */
494            private final ValueAndJacobianFunction model;
495            /** Target values for the model function at optimum. */
496            private final RealVector target;
497
498            /**
499             * Create an {@link Evaluation} with no weights.
500             *
501             * @param model  the model function
502             * @param target the observed values
503             * @param point  the abscissa
504             */
505            private LazyUnweightedEvaluation(final ValueAndJacobianFunction model,
506                                             final RealVector target,
507                                             final RealVector point) {
508                super(target.getDimension());
509                // Safe to cast as long as we control usage of this class.
510                this.model = model;
511                this.point = point;
512                this.target = target;
513            }
514
515            /** {@inheritDoc} */
516            public RealMatrix getJacobian() {
517                return model.computeJacobian(point.toArray());
518            }
519
520            /** {@inheritDoc} */
521            public RealVector getPoint() {
522                return point;
523            }
524
525            /** {@inheritDoc} */
526            public RealVector getResiduals() {
527                return target.subtract(model.computeValue(point.toArray()));
528            }
529        }
530    }
531}
532