UnscentedTransformProvider.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.util;

  18. import org.hipparchus.linear.ArrayRealVector;
  19. import org.hipparchus.linear.MatrixUtils;
  20. import org.hipparchus.linear.RealMatrix;
  21. import org.hipparchus.linear.RealVector;

  22. /**
  23.  * Provider for unscented transform.
  24.  * @since 2.2
  25.  */
  26. public interface UnscentedTransformProvider {

  27.     /**
  28.      * Perform the unscented transform from a state and its covariance.
  29.      * @param state process state
  30.      * @param covariance covariance associated with the process state
  31.      * @return an array containing the sigma points of the unscented transform
  32.      */
  33.     RealVector[] unscentedTransform(RealVector state, RealMatrix covariance);

  34.     /**
  35.      * Computes a weighted mean state from a given set of sigma points.
  36.      * <p>
  37.      * This method can be used for computing both the mean state and the mean measurement
  38.      * in an Unscented Kalman filter.
  39.      * </p>
  40.      * <p>
  41.      * It corresponds to Equation 17 of "Wan, E. A., &amp; Van Der Merwe, R. The unscented Kalman filter for nonlinear estimation"
  42.      * </p>
  43.      * @param sigmaPoints input samples
  44.      * @return weighted mean state
  45.      */
  46.     default RealVector getUnscentedMeanState(RealVector[] sigmaPoints) {

  47.         // Sigma point dimension
  48.         final int sigmaPointDimension = sigmaPoints[0].getDimension();

  49.         // Compute weighted mean
  50.         // ---------------------

  51.         RealVector weightedMean = new ArrayRealVector(sigmaPointDimension);

  52.         // Compute the weight coefficients wm
  53.         final RealVector wm = getWm();

  54.         // Weight each sigma point and sum them
  55.         for (int i = 0; i < sigmaPoints.length; i++) {
  56.             weightedMean = weightedMean.add(sigmaPoints[i].mapMultiply(wm.getEntry(i)));
  57.         }

  58.         return weightedMean;
  59.     }

  60.     /** Computes the unscented covariance matrix from a weighted mean state and a set of sigma points.
  61.      * <p>
  62.      * This method can be used for computing both the predicted state
  63.      * covariance matrix and the innovation covariance matrix in an Unscented Kalman filter.
  64.      * </p>
  65.      * <p>
  66.      * It corresponds to Equation 18 of "Wan, E. A., &amp; Van Der Merwe, R. The unscented Kalman filter for nonlinear estimation"
  67.      * </p>
  68.      * @param sigmaPoints input sigma points
  69.      * @param meanState weighted mean state
  70.      * @return the unscented covariance matrix
  71.      */
  72.     default RealMatrix getUnscentedCovariance(RealVector[] sigmaPoints, RealVector meanState) {

  73.         // State dimension
  74.         final int stateDimension = meanState.getDimension();

  75.         // Compute covariance matrix
  76.         // -------------------------

  77.         RealMatrix covarianceMatrix = MatrixUtils.createRealMatrix(stateDimension, stateDimension);

  78.         // Compute the weight coefficients wc
  79.         final RealVector wc = getWc();

  80.         // Reconstruct the covariance
  81.         for (int i = 0; i < sigmaPoints.length; i++) {
  82.             final RealMatrix diff = MatrixUtils.createColumnRealMatrix(sigmaPoints[i].subtract(meanState).toArray());
  83.             covarianceMatrix = covarianceMatrix.add(diff.multiplyTransposed(diff).scalarMultiply(wc.getEntry(i)));
  84.         }

  85.         return covarianceMatrix;
  86.     }

  87.     /**
  88.      * Perform the inverse unscented transform from an array of sigma points.
  89.      * @param sigmaPoints array containing the sigma points of the unscented transform
  90.      * @return mean state and associated covariance
  91.      */
  92.     default Pair<RealVector, RealMatrix> inverseUnscentedTransform(RealVector[] sigmaPoints) {

  93.         // Mean state
  94.         final RealVector meanState = getUnscentedMeanState(sigmaPoints);

  95.         // Return state and covariance
  96.         return new Pair<>(meanState, getUnscentedCovariance(sigmaPoints, meanState));
  97.     }

  98.     /**
  99.      * Get the covariance weights.
  100.      * @return the covariance weights
  101.      */
  102.     RealVector getWc();

  103.     /**
  104.      * Get the mean weights.
  105.      * @return the mean weights
  106.      */
  107.     RealVector getWm();

  108. }