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 /*
19 * This is not the original file distributed by the Apache Software Foundation
20 * It has been modified by the Hipparchus project
21 */
22 package org.hipparchus.optim.nonlinear.vector.leastsquares;
23
24 import org.hipparchus.linear.ArrayRealVector;
25 import org.hipparchus.linear.DecompositionSolver;
26 import org.hipparchus.linear.QRDecomposition;
27 import org.hipparchus.linear.RealMatrix;
28 import org.hipparchus.linear.RealVector;
29 import org.hipparchus.optim.nonlinear.vector.leastsquares.LeastSquaresProblem.Evaluation;
30 import org.hipparchus.util.FastMath;
31
32 /**
33 * An implementation of {@link Evaluation} that is designed for extension. All of the
34 * methods implemented here use the methods that are left unimplemented.
35 */
36 public abstract class AbstractEvaluation implements Evaluation {
37
38 /** number of observations */
39 private final int observationSize;
40
41 /**
42 * Constructor.
43 *
44 * @param observationSize the number of observations.
45 * Needed for {@link #getRMS()} and {@link #getReducedChiSquare(int)}.
46 */
47 protected AbstractEvaluation(final int observationSize) {
48 this.observationSize = observationSize;
49 }
50
51 /** {@inheritDoc} */
52 @Override
53 public RealMatrix getCovariances(double threshold) {
54 // Set up the Jacobian.
55 final RealMatrix j = this.getJacobian();
56
57 // Compute transpose(J)J.
58 final RealMatrix jTj = j.transposeMultiply(j);
59
60 // Compute the covariances matrix.
61 final DecompositionSolver solver
62 = new QRDecomposition(jTj, threshold).getSolver();
63 return solver.getInverse();
64 }
65
66 /** {@inheritDoc} */
67 @Override
68 public RealVector getSigma(double covarianceSingularityThreshold) {
69 final RealMatrix cov = this.getCovariances(covarianceSingularityThreshold);
70 final int nC = cov.getColumnDimension();
71 final RealVector sig = new ArrayRealVector(nC);
72 for (int i = 0; i < nC; ++i) {
73 sig.setEntry(i, FastMath.sqrt(cov.getEntry(i,i)));
74 }
75 return sig;
76 }
77
78 /** {@inheritDoc} */
79 @Override
80 public double getRMS() {
81 return FastMath.sqrt(getReducedChiSquare(1));
82 }
83
84 /** {@inheritDoc} */
85 @Override
86 public double getCost() {
87 return FastMath.sqrt(getChiSquare());
88 }
89
90 /** {@inheritDoc} */
91 @Override
92 public double getChiSquare() {
93 final ArrayRealVector r = new ArrayRealVector(getResiduals());
94 return r.dotProduct(r);
95 }
96
97 /** {@inheritDoc} */
98 @Override
99 public double getReducedChiSquare(int numberOfFittedParameters) {
100 return getChiSquare() / (observationSize - numberOfFittedParameters + 1);
101 }
102 }