JacobiansMapper.java

  1. /* Copyright 2002-2019 CS Systèmes d'Information
  2.  * Licensed to CS Systèmes d'Information (CS) under one or more
  3.  * contributor license agreements.  See the NOTICE file distributed with
  4.  * this work for additional information regarding copyright ownership.
  5.  * CS 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.  *   http://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.orekit.propagation.numerical;

  18. import org.hipparchus.linear.Array2DRowRealMatrix;
  19. import org.hipparchus.linear.DecompositionSolver;
  20. import org.hipparchus.linear.QRDecomposition;
  21. import org.hipparchus.linear.RealMatrix;
  22. import org.orekit.orbits.Orbit;
  23. import org.orekit.orbits.OrbitType;
  24. import org.orekit.orbits.PositionAngle;
  25. import org.orekit.propagation.SpacecraftState;
  26. import org.orekit.utils.ParameterDriversList;

  27. /** Mapper between two-dimensional Jacobian matrices and one-dimensional {@link
  28.  * SpacecraftState#getAdditionalState(String) additional state arrays}.
  29.  * <p>
  30.  * This class does not hold the states by itself. Instances of this class are guaranteed
  31.  * to be immutable.
  32.  * </p>
  33.  * @author Luc Maisonobe
  34.  * @see org.orekit.propagation.numerical.PartialDerivativesEquations
  35.  * @see org.orekit.propagation.numerical.NumericalPropagator
  36.  * @see SpacecraftState#getAdditionalState(String)
  37.  * @see org.orekit.propagation.AbstractPropagator
  38.  */
  39. public class JacobiansMapper {

  40.     /** State dimension, fixed to 6.
  41.      * @since 9.0
  42.      */
  43.     public static final int STATE_DIMENSION = 6;

  44.     /** Name. */
  45.     private String name;

  46.     /** Selected parameters for Jacobian computation. */
  47.     private final ParameterDriversList parameters;

  48.     /** Orbit type. */
  49.     private final OrbitType orbitType;

  50.     /** Position angle type. */
  51.     private final PositionAngle angleType;

  52.     /** Simple constructor.
  53.      * @param name name of the Jacobians
  54.      * @param parameters selected parameters for Jacobian computation
  55.      * @param orbitType orbit type
  56.      * @param angleType position angle type
  57.      */
  58.     JacobiansMapper(final String name, final ParameterDriversList parameters,
  59.                     final OrbitType orbitType, final PositionAngle angleType) {
  60.         this.name       = name;
  61.         this.parameters = parameters;
  62.         this.orbitType  = orbitType;
  63.         this.angleType  = angleType;
  64.     }

  65.     /** Get the name of the partial Jacobians.
  66.      * @return name of the Jacobians
  67.      */
  68.     public String getName() {
  69.         return name;
  70.     }

  71.     /** Compute the length of the one-dimensional additional state array needed.
  72.      * @return length of the one-dimensional additional state array
  73.      */
  74.     public int getAdditionalStateDimension() {
  75.         return STATE_DIMENSION * (STATE_DIMENSION + parameters.getNbParams());
  76.     }

  77.     /** Get the state vector dimension.
  78.      * @return state vector dimension
  79.      * @deprecated as of 9.0, replaced with {@link #STATE_DIMENSION}
  80.      */
  81.     @Deprecated
  82.     public int getStateDimension() {
  83.         return STATE_DIMENSION;
  84.     }

  85.     /** Get the number of parameters.
  86.      * @return number of parameters
  87.      */
  88.     public int getParameters() {
  89.         return parameters.getNbParams();
  90.     }

  91.     /** Get the conversion Jacobian between state parameters and Cartesian parameters.
  92.      * @param state spacecraft state
  93.      * @return conversion Jacobian
  94.      */
  95.     private double[][] getdYdC(final SpacecraftState state) {

  96.         final double[][] dYdC = new double[STATE_DIMENSION][STATE_DIMENSION];

  97.         // make sure the state is in the desired orbit type
  98.         final Orbit orbit = orbitType.convertType(state.getOrbit());

  99.         // compute the Jacobian, taking the position angle type into account
  100.         orbit.getJacobianWrtCartesian(angleType, dYdC);

  101.         return dYdC;

  102.     }

  103.     /** Set the Jacobian with respect to state into a one-dimensional additional state array.
  104.      * <p>
  105.      * This method converts the Jacobians to Cartesian parameters and put the converted data
  106.      * in the one-dimensional {@code p} array.
  107.      * </p>
  108.      * @param state spacecraft state
  109.      * @param dY1dY0 Jacobian of current state at time t₁
  110.      * with respect to state at some previous time t₀
  111.      * @param dY1dP Jacobian of current state at time t₁
  112.      * with respect to parameters (may be null if there are no parameters)
  113.      * @param p placeholder where to put the one-dimensional additional state
  114.      * @see #getStateJacobian(SpacecraftState, double[][])
  115.      */
  116.     void setInitialJacobians(final SpacecraftState state, final double[][] dY1dY0,
  117.                              final double[][] dY1dP, final double[] p) {

  118.         // set up a converter between state parameters and Cartesian parameters
  119.         final RealMatrix dY1dC1 = new Array2DRowRealMatrix(getdYdC(state), false);
  120.         final DecompositionSolver solver = new QRDecomposition(dY1dC1).getSolver();

  121.         // convert the provided state Jacobian to Cartesian parameters
  122.         final RealMatrix dC1dY0 = solver.solve(new Array2DRowRealMatrix(dY1dY0, false));

  123.         // map the converted state Jacobian to one-dimensional array
  124.         int index = 0;
  125.         for (int i = 0; i < STATE_DIMENSION; ++i) {
  126.             for (int j = 0; j < STATE_DIMENSION; ++j) {
  127.                 p[index++] = dC1dY0.getEntry(i, j);
  128.             }
  129.         }

  130.         if (parameters.getNbParams() != 0) {
  131.             // convert the provided state Jacobian to Cartesian parameters
  132.             final RealMatrix dC1dP = solver.solve(new Array2DRowRealMatrix(dY1dP, false));

  133.             // map the converted parameters Jacobian to one-dimensional array
  134.             for (int i = 0; i < STATE_DIMENSION; ++i) {
  135.                 for (int j = 0; j < parameters.getNbParams(); ++j) {
  136.                     p[index++] = dC1dP.getEntry(i, j);
  137.                 }
  138.             }
  139.         }

  140.     }

  141.     /** Get the Jacobian with respect to state from a one-dimensional additional state array.
  142.      * <p>
  143.      * This method extract the data from the {@code state} and put it in the
  144.      * {@code dYdY0} array.
  145.      * </p>
  146.      * @param state spacecraft state
  147.      * @param dYdY0 placeholder where to put the Jacobian with respect to state
  148.           * @see #getParametersJacobian(SpacecraftState, double[][])
  149.      */
  150.     public void getStateJacobian(final SpacecraftState state,  final double[][] dYdY0) {

  151.         // get the conversion Jacobian between state parameters and Cartesian parameters
  152.         final double[][] dYdC = getdYdC(state);

  153.         // extract the additional state
  154.         final double[] p = state.getAdditionalState(name);

  155.         // compute dYdY0 = dYdC * dCdY0, without allocating new arrays
  156.         for (int i = 0; i < STATE_DIMENSION; i++) {
  157.             final double[] rowC = dYdC[i];
  158.             final double[] rowD = dYdY0[i];
  159.             for (int j = 0; j < STATE_DIMENSION; ++j) {
  160.                 double sum = 0;
  161.                 int pIndex = j;
  162.                 for (int k = 0; k < STATE_DIMENSION; ++k) {
  163.                     sum += rowC[k] * p[pIndex];
  164.                     pIndex += STATE_DIMENSION;
  165.                 }
  166.                 rowD[j] = sum;
  167.             }
  168.         }

  169.     }

  170.     /** Get theJacobian with respect to parameters from a one-dimensional additional state array.
  171.      * <p>
  172.      * This method extract the data from the {@code state} and put it in the
  173.      * {@code dYdP} array.
  174.      * </p>
  175.      * <p>
  176.      * If no parameters have been set in the constructor, the method returns immediately and
  177.      * does not reference {@code dYdP} which can safely be null in this case.
  178.      * </p>
  179.      * @param state spacecraft state
  180.      * @param dYdP placeholder where to put the Jacobian with respect to parameters
  181.           * @see #getStateJacobian(SpacecraftState, double[][])
  182.      */
  183.     public void getParametersJacobian(final SpacecraftState state, final double[][] dYdP) {

  184.         if (parameters.getNbParams() != 0) {

  185.             // get the conversion Jacobian between state parameters and Cartesian parameters
  186.             final double[][] dYdC = getdYdC(state);

  187.             // extract the additional state
  188.             final double[] p = state.getAdditionalState(name);

  189.             // compute dYdP = dYdC * dCdP, without allocating new arrays
  190.             for (int i = 0; i < STATE_DIMENSION; i++) {
  191.                 final double[] rowC = dYdC[i];
  192.                 final double[] rowD = dYdP[i];
  193.                 for (int j = 0; j < parameters.getNbParams(); ++j) {
  194.                     double sum = 0;
  195.                     int pIndex = j + STATE_DIMENSION * STATE_DIMENSION;
  196.                     for (int k = 0; k < STATE_DIMENSION; ++k) {
  197.                         sum += rowC[k] * p[pIndex];
  198.                         pIndex += parameters.getNbParams();
  199.                     }
  200.                     rowD[j] = sum;
  201.                 }
  202.             }

  203.         }

  204.     }

  205. }