AbstractKalmanFilter.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.filtering.kalman;

  18. import org.hipparchus.exception.MathIllegalArgumentException;
  19. import org.hipparchus.linear.MatrixDecomposer;
  20. import org.hipparchus.linear.RealMatrix;
  21. import org.hipparchus.linear.RealVector;

  22. /**
  23.  * Shared parts between linear and non-linear Kalman filters.
  24.  * @param <T> the type of the measurements
  25.  * @since 1.3
  26.  */
  27. public abstract class AbstractKalmanFilter<T extends Measurement> implements KalmanFilter<T> {

  28.     /** Decomposer decomposer to use for the correction phase. */
  29.     private final MatrixDecomposer decomposer;

  30.     /** Predicted state. */
  31.     private ProcessEstimate predicted;

  32.     /** Corrected state. */
  33.     private ProcessEstimate corrected;

  34.     /** Prior corrected covariance. */
  35.     private RealMatrix priorCovariance;

  36.     /** State transition matrix. */
  37.     private RealMatrix stateTransitionMatrix;

  38.     /** Observer. */
  39.     private KalmanObserver observer;


  40.     /** Simple constructor.
  41.      * @param decomposer decomposer to use for the correction phase
  42.      * @param initialState initial state
  43.      */
  44.     protected AbstractKalmanFilter(final MatrixDecomposer decomposer, final ProcessEstimate initialState) {
  45.         this.decomposer = decomposer;
  46.         this.corrected  = initialState;
  47.         this.priorCovariance = null;
  48.         this.stateTransitionMatrix = null;
  49.         this.observer = null;
  50.     }

  51.     /** Perform prediction step.
  52.      * @param time process time
  53.      * @param predictedState predicted state vector
  54.      * @param stm state transition matrix
  55.      * @param noise process noise covariance matrix
  56.      */
  57.     protected void predict(final double time, final RealVector predictedState, final RealMatrix stm, final RealMatrix noise) {
  58.         final RealMatrix predictedCovariance =
  59.                         stm.multiply(corrected.getCovariance().multiplyTransposed(stm)).add(noise);
  60.         predicted = new ProcessEstimate(time, predictedState, predictedCovariance);
  61.         stateTransitionMatrix = stm;
  62.         priorCovariance = corrected.getCovariance();
  63.         corrected = null;
  64.     }

  65.     /** Compute innovation covariance matrix.
  66.      * @param r measurement covariance
  67.      * @param h Jacobian of the measurement with respect to the state
  68.      * (may be null if measurement should be ignored)
  69.      * @return innovation covariance matrix, defined as \(h.P.h^T + r\), or
  70.      * null if h is null
  71.      */
  72.     protected RealMatrix computeInnovationCovarianceMatrix(final RealMatrix r, final RealMatrix h) {
  73.         if (h == null) {
  74.             return null;
  75.         }
  76.         final RealMatrix phT = predicted.getCovariance().multiplyTransposed(h);
  77.         return h.multiply(phT).add(r);
  78.     }

  79.     /** Perform correction step.
  80.      * @param measurement single measurement to handle
  81.      * @param stm state transition matrix
  82.      * @param innovation innovation vector (i.e. residuals)
  83.      * (may be null if measurement should be ignored)
  84.      * @param h Jacobian of the measurement with respect to the state
  85.      * (may be null if measurement should be ignored)
  86.      * @param s innovation covariance matrix
  87.      * (may be null if measurement should be ignored)
  88.      * @exception MathIllegalArgumentException if matrix cannot be decomposed
  89.      */
  90.     protected void correct(final T measurement, final RealMatrix stm, final RealVector innovation,
  91.                            final RealMatrix h, final RealMatrix s)
  92.         throws MathIllegalArgumentException {

  93.         if (innovation == null) {
  94.             // measurement should be ignored
  95.             corrected = predicted;
  96.             return;
  97.         }

  98.         // compute Kalman gain k
  99.         // the following is equivalent to k = p.h^T * (h.p.h^T + r)^(-1)
  100.         // we don't want to compute the inverse of a matrix,
  101.         // we start by post-multiplying by h.p.h^T + r and get
  102.         // k.(h.p.h^T + r) = p.h^T
  103.         // then we transpose, knowing that both p and r are symmetric matrices
  104.         // (h.p.h^T + r).k^T = h.p
  105.         // then we can use linear system solving instead of matrix inversion
  106.         final RealMatrix k = decomposer.
  107.                              decompose(s).
  108.                              solve(h.multiply(predicted.getCovariance())).
  109.                              transpose();

  110.         // correct state vector
  111.         final RealVector correctedState = predicted.getState().add(k.operate(innovation));

  112.         // here we use the Joseph algorithm (see "Fundamentals of Astrodynamics and Applications,
  113.         // Vallado, Fourth Edition ยง10.6 eq.10-34) which is equivalent to
  114.         // the traditional Pest = (I - k.h) x Ppred expression but guarantees the output stays symmetric:
  115.         // Pest = (I -k.h) Ppred (I - k.h)^T + k.r.k^T
  116.         final RealMatrix idMkh = k.multiply(h);
  117.         for (int i = 0; i < idMkh.getRowDimension(); ++i) {
  118.             for (int j = 0; j < idMkh.getColumnDimension(); ++j) {
  119.                 idMkh.multiplyEntry(i, j, -1);
  120.             }
  121.             idMkh.addToEntry(i, i, 1.0);
  122.         }
  123.         final RealMatrix r = measurement.getCovariance();
  124.         final RealMatrix correctedCovariance =
  125.                         idMkh.multiply(predicted.getCovariance()).multiplyTransposed(idMkh).
  126.                         add(k.multiply(r).multiplyTransposed(k));

  127.         corrected = new ProcessEstimate(measurement.getTime(), correctedState, correctedCovariance,
  128.                                         stm, h, s, k);

  129.     }

  130.     /** Get the observer.
  131.      * @return the observer
  132.      */
  133.     protected KalmanObserver getObserver() {
  134.         return observer;
  135.     }


  136.     /** {@inheritDoc} */
  137.     @Override
  138.     public void setObserver(final KalmanObserver kalmanObserver) {
  139.         observer = kalmanObserver;
  140.         observer.init(this);
  141.     }

  142.     /** Get the predicted state.
  143.      * @return predicted state
  144.      */
  145.     @Override
  146.     public ProcessEstimate getPredicted() {
  147.         return predicted;
  148.     }

  149.     /** Get the corrected state.
  150.      * @return corrected state
  151.      */
  152.     @Override
  153.     public ProcessEstimate getCorrected() {
  154.         return corrected;
  155.     }

  156.     /** {@inheritDoc} */
  157.     @Override
  158.     public RealMatrix getStateCrossCovariance() {
  159.         return priorCovariance.multiplyTransposed(stateTransitionMatrix);
  160.     }
  161. }