EpochDerivativesEquations.java

  1. /* Copyright 2002-2022 CS GROUP
  2.  * Licensed to CS GROUP (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 java.util.IdentityHashMap;
  19. import java.util.Map;

  20. import org.hipparchus.analysis.differentiation.Gradient;
  21. import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
  22. import org.orekit.errors.OrekitException;
  23. import org.orekit.errors.OrekitMessages;
  24. import org.orekit.forces.ForceModel;
  25. import org.orekit.forces.gravity.ThirdBodyAttractionEpoch;
  26. import org.orekit.propagation.FieldSpacecraftState;
  27. import org.orekit.propagation.SpacecraftState;
  28. import org.orekit.propagation.integration.AdditionalDerivativesProvider;
  29. import org.orekit.time.AbsoluteDate;
  30. import org.orekit.utils.ParameterDriver;
  31. import org.orekit.utils.ParameterDriversList;

  32. /** This class is a copy of {@link AbsolutePartialDerivativesEquations}
  33.  *  The computation of the derivatives of the acceleration due to a ThirdBodyAttraction
  34.  *  has been added.
  35.  *
  36.  * {@link AdditionalDerivativesProvider Provider} computing the partial derivatives
  37.  * of the state (orbit) with respect to initial state and force models parameters.
  38.  * <p>
  39.  * This set of equations are automatically added to a {@link NumericalPropagator numerical propagator}
  40.  * in order to compute partial derivatives of the orbit along with the orbit itself. This is
  41.  * useful for example in orbit determination applications.
  42.  * </p>
  43.  * <p>
  44.  * The partial derivatives with respect to initial state can be either dimension 6
  45.  * (orbit only) or 7 (orbit and mass).
  46.  * </p>
  47.  * <p>
  48.  * The partial derivatives with respect to force models parameters has a dimension
  49.  * equal to the number of selected parameters. Parameters selection is implemented at
  50.  * {@link ForceModel force models} level. Users must retrieve a {@link ParameterDriver
  51.  * parameter driver} using {@link ForceModel#getParameterDriver(String)} and then
  52.  * select it by calling {@link ParameterDriver#setSelected(boolean) setSelected(true)}.
  53.  * </p>
  54.  * <p>
  55.  * If several force models provide different {@link ParameterDriver drivers} for the
  56.  * same parameter name, selecting any of these drivers has the side effect of
  57.  * selecting all the drivers for this shared parameter. In this case, the partial
  58.  * derivatives will be the sum of the partial derivatives contributed by the
  59.  * corresponding force models. This case typically arises for central attraction
  60.  * coefficient, which has an influence on {@link org.orekit.forces.gravity.NewtonianAttraction
  61.  * Newtonian attraction}, {@link org.orekit.forces.gravity.HolmesFeatherstoneAttractionModel
  62.  * gravity field}, and {@link org.orekit.forces.gravity.Relativity relativity}.
  63.  * </p>
  64.  * @author V&eacute;ronique Pommier-Maurussane
  65.  * @author Luc Maisonobe
  66.  * @since 10.2
  67.  */
  68. @SuppressWarnings("deprecation")
  69. public class EpochDerivativesEquations
  70.     implements AdditionalDerivativesProvider,
  71.     org.orekit.propagation.integration.AdditionalEquations  {

  72.     /** Propagator computing state evolution. */
  73.     private final NumericalPropagator propagator;

  74.     /** Selected parameters for Jacobian computation. */
  75.     private ParameterDriversList selected;

  76.     /** Parameters map. */
  77.     private Map<ParameterDriver, Integer> map;

  78.     /** Name. */
  79.     private final String name;

  80.     /** Flag for Jacobian matrices initialization. */
  81.     private boolean initialized;

  82.     /** Simple constructor.
  83.      * <p>
  84.      * Upon construction, this set of equations is <em>automatically</em> added to
  85.      * the propagator by calling its {@link
  86.      * NumericalPropagator#addAdditionalEquations(AdditionalEquations)} method. So
  87.      * there is no need to call this method explicitly for these equations.
  88.      * </p>
  89.      * @param name name of the partial derivatives equations
  90.      * @param propagator the propagator that will handle the orbit propagation
  91.      */
  92.     public EpochDerivativesEquations(final String name, final NumericalPropagator propagator) {
  93.         this.name                   = name;
  94.         this.selected               = null;
  95.         this.map                    = null;
  96.         this.propagator             = propagator;
  97.         this.initialized            = false;
  98.         propagator.addAdditionalDerivativesProvider(this);
  99.     }

  100.     /** {@inheritDoc} */
  101.     public String getName() {
  102.         return name;
  103.     }

  104.     /** {@inheritDoc} */
  105.     @Override
  106.     public int getDimension() {
  107.         freezeParametersSelection();
  108.         return 6 * (6 + selected.getNbParams() + 1);
  109.     }

  110.     /** Freeze the selected parameters from the force models.
  111.      */
  112.     private void freezeParametersSelection() {
  113.         if (selected == null) {

  114.             // first pass: gather all parameters, binding similar names together
  115.             selected = new ParameterDriversList();
  116.             for (final ForceModel provider : propagator.getAllForceModels()) {
  117.                 for (final ParameterDriver driver : provider.getParametersDrivers()) {
  118.                     selected.add(driver);
  119.                 }
  120.             }

  121.             // second pass: now that shared parameter names are bound together,
  122.             // their selections status have been synchronized, we can filter them
  123.             selected.filter(true);

  124.             // third pass: sort parameters lexicographically
  125.             selected.sort();

  126.             // fourth pass: set up a map between parameters drivers and matrices columns
  127.             map = new IdentityHashMap<>();
  128.             int parameterIndex = 0;
  129.             for (final ParameterDriver selectedDriver : selected.getDrivers()) {
  130.                 for (final ForceModel provider : propagator.getAllForceModels()) {
  131.                     for (final ParameterDriver driver : provider.getParametersDrivers()) {
  132.                         if (driver.getName().equals(selectedDriver.getName())) {
  133.                             map.put(driver, parameterIndex);
  134.                         }
  135.                     }
  136.                 }
  137.                 ++parameterIndex;
  138.             }

  139.         }
  140.     }

  141.     /** Get the selected parameters, in Jacobian matrix column order.
  142.      * <p>
  143.      * The force models parameters for which partial derivatives are desired,
  144.      * <em>must</em> have been {@link ParameterDriver#setSelected(boolean) selected}
  145.      * before this method is called, so the proper list is returned.
  146.      * </p>
  147.      * @return selected parameters, in Jacobian matrix column order which
  148.      * is lexicographic order
  149.      */
  150.     public ParameterDriversList getSelectedParameters() {
  151.         freezeParametersSelection();
  152.         return selected;
  153.     }

  154.     /** Set the initial value of the Jacobian with respect to state and parameter.
  155.      * <p>
  156.      * This method is equivalent to call {@link #setInitialJacobians(SpacecraftState,
  157.      * double[][], double[][])} with dYdY0 set to the identity matrix and dYdP set
  158.      * to a zero matrix.
  159.      * </p>
  160.      * <p>
  161.      * The force models parameters for which partial derivatives are desired,
  162.      * <em>must</em> have been {@link ParameterDriver#setSelected(boolean) selected}
  163.      * before this method is called, so proper matrices dimensions are used.
  164.      * </p>
  165.      * @param s0 initial state
  166.      * @return state with initial Jacobians added
  167.      * @see #getSelectedParameters()
  168.      * @since 9.0
  169.      */
  170.     public SpacecraftState setInitialJacobians(final SpacecraftState s0) {
  171.         freezeParametersSelection();
  172.         final int epochStateDimension = 6;
  173.         final double[][] dYdY0 = new double[epochStateDimension][epochStateDimension];
  174.         final double[][] dYdP  = new double[epochStateDimension][selected.getNbParams() + 6];
  175.         for (int i = 0; i < epochStateDimension; ++i) {
  176.             dYdY0[i][i] = 1.0;
  177.         }
  178.         return setInitialJacobians(s0, dYdY0, dYdP);
  179.     }

  180.     /** Set the initial value of the Jacobian with respect to state and parameter.
  181.      * <p>
  182.      * The returned state must be added to the propagator (it is not done
  183.      * automatically, as the user may need to add more states to it).
  184.      * </p>
  185.      * <p>
  186.      * The force models parameters for which partial derivatives are desired,
  187.      * <em>must</em> have been {@link ParameterDriver#setSelected(boolean) selected}
  188.      * before this method is called, and the {@code dY1dP} matrix dimension <em>must</em>
  189.      * be consistent with the selection.
  190.      * </p>
  191.      * @param s1 current state
  192.      * @param dY1dY0 Jacobian of current state at time t₁ with respect
  193.      * to state at some previous time t₀ (must be 6x6)
  194.      * @param dY1dP Jacobian of current state at time t₁ with respect
  195.      * to parameters (may be null if no parameters are selected)
  196.      * @return state with initial Jacobians added
  197.      * @see #getSelectedParameters()
  198.      */
  199.     public SpacecraftState setInitialJacobians(final SpacecraftState s1,
  200.                                                final double[][] dY1dY0, final double[][] dY1dP) {

  201.         freezeParametersSelection();

  202.         // Check dimensions
  203.         final int stateDimEpoch = dY1dY0.length;
  204.         if (stateDimEpoch != 6 || stateDimEpoch != dY1dY0[0].length) {
  205.             throw new OrekitException(OrekitMessages.STATE_JACOBIAN_NOT_6X6,
  206.                                       stateDimEpoch, dY1dY0[0].length);
  207.         }
  208.         if (dY1dP != null && stateDimEpoch != dY1dP.length) {
  209.             throw new OrekitException(OrekitMessages.STATE_AND_PARAMETERS_JACOBIANS_ROWS_MISMATCH,
  210.                                       stateDimEpoch, dY1dP.length);
  211.         }

  212.         // store the matrices as a single dimension array
  213.         initialized = true;
  214.         final AbsoluteJacobiansMapper absoluteMapper = getMapper();
  215.         final double[] p = new double[absoluteMapper.getAdditionalStateDimension() + 6];
  216.         absoluteMapper.setInitialJacobians(s1, dY1dY0, dY1dP, p);

  217.         // set value in propagator
  218.         return s1.addAdditionalState(name, p);

  219.     }

  220.     /** Get a mapper between two-dimensional Jacobians and one-dimensional additional state.
  221.      * @return a mapper between two-dimensional Jacobians and one-dimensional additional state,
  222.      * with the same name as the instance
  223.      * @see #setInitialJacobians(SpacecraftState)
  224.      * @see #setInitialJacobians(SpacecraftState, double[][], double[][])
  225.      */
  226.     public AbsoluteJacobiansMapper getMapper() {
  227.         if (!initialized) {
  228.             throw new OrekitException(OrekitMessages.STATE_JACOBIAN_NOT_INITIALIZED);
  229.         }
  230.         return new AbsoluteJacobiansMapper(name, selected);
  231.     }

  232.     /** {@inheritDoc} */
  233.     public void init(final SpacecraftState initialState, final AbsoluteDate target) {
  234.         // FIXME: remove in 12.0 when AdditionalEquations is removed
  235.         AdditionalDerivativesProvider.super.init(initialState, target);
  236.     }

  237.     /** {@inheritDoc} */
  238.     public double[] computeDerivatives(final SpacecraftState s, final double[] pDot) {
  239.         // FIXME: remove in 12.0 when AdditionalEquations is removed
  240.         System.arraycopy(derivatives(s), 0, pDot, 0, pDot.length);
  241.         return null;
  242.     }

  243.     /** {@inheritDoc} */
  244.     public double[] derivatives(final SpacecraftState s) {

  245.         // initialize acceleration Jacobians to zero
  246.         final int paramDimEpoch = selected.getNbParams() + 1; // added epoch
  247.         final int dimEpoch      = 3;
  248.         final double[][] dAccdParam = new double[dimEpoch][paramDimEpoch];
  249.         final double[][] dAccdPos   = new double[dimEpoch][dimEpoch];
  250.         final double[][] dAccdVel   = new double[dimEpoch][dimEpoch];

  251.         final NumericalGradientConverter fullConverter    = new NumericalGradientConverter(s, 6, propagator.getAttitudeProvider());
  252.         final NumericalGradientConverter posOnlyConverter = new NumericalGradientConverter(s, 3, propagator.getAttitudeProvider());

  253.         // compute acceleration Jacobians, finishing with the largest force: Newtonian attraction
  254.         for (final ForceModel forceModel : propagator.getAllForceModels()) {
  255.             final NumericalGradientConverter converter = forceModel.dependsOnPositionOnly() ? posOnlyConverter : fullConverter;
  256.             final FieldSpacecraftState<Gradient> dsState = converter.getState(forceModel);
  257.             final Gradient[] parameters = converter.getParameters(dsState, forceModel);

  258.             final FieldVector3D<Gradient> acceleration = forceModel.acceleration(dsState, parameters);
  259.             final double[] derivativesX = acceleration.getX().getGradient();
  260.             final double[] derivativesY = acceleration.getY().getGradient();
  261.             final double[] derivativesZ = acceleration.getZ().getGradient();

  262.             // update Jacobians with respect to state
  263.             addToRow(derivativesX, 0, converter.getFreeStateParameters(), dAccdPos, dAccdVel);
  264.             addToRow(derivativesY, 1, converter.getFreeStateParameters(), dAccdPos, dAccdVel);
  265.             addToRow(derivativesZ, 2, converter.getFreeStateParameters(), dAccdPos, dAccdVel);

  266.             int index = converter.getFreeStateParameters();
  267.             for (ParameterDriver driver : forceModel.getParametersDrivers()) {
  268.                 if (driver.isSelected()) {
  269.                     final int parameterIndex = map.get(driver);
  270.                     dAccdParam[0][parameterIndex] += derivativesX[index];
  271.                     dAccdParam[1][parameterIndex] += derivativesY[index];
  272.                     dAccdParam[2][parameterIndex] += derivativesZ[index];
  273.                     ++index;
  274.                 }
  275.             }

  276.             // Add the derivatives of the acceleration w.r.t. the Epoch
  277.             if (forceModel instanceof ThirdBodyAttractionEpoch) {
  278.                 final double[] parametersValues = new double[] {parameters[0].getValue()};
  279.                 final double[] derivatives = ((ThirdBodyAttractionEpoch) forceModel).getDerivativesToEpoch(s, parametersValues);
  280.                 dAccdParam[0][paramDimEpoch - 1] += derivatives[0];
  281.                 dAccdParam[1][paramDimEpoch - 1] += derivatives[1];
  282.                 dAccdParam[2][paramDimEpoch - 1] += derivatives[2];
  283.             }

  284.         }

  285.         // the variational equations of the complete state Jacobian matrix have the following form:

  286.         // [        |        ]   [                 |                  ]   [     |     ]
  287.         // [  Adot  |  Bdot  ]   [  dVel/dPos = 0  |  dVel/dVel = Id  ]   [  A  |  B  ]
  288.         // [        |        ]   [                 |                  ]   [     |     ]
  289.         // ---------+---------   ------------------+------------------- * ------+------
  290.         // [        |        ]   [                 |                  ]   [     |     ]
  291.         // [  Cdot  |  Ddot  ] = [    dAcc/dPos    |     dAcc/dVel    ]   [  C  |  D  ]
  292.         // [        |        ]   [                 |                  ]   [     |     ]

  293.         // The A, B, C and D sub-matrices and their derivatives (Adot ...) are 3x3 matrices

  294.         // The expanded multiplication above can be rewritten to take into account
  295.         // the fixed values found in the sub-matrices in the left factor. This leads to:

  296.         //     [ Adot ] = [ C ]
  297.         //     [ Bdot ] = [ D ]
  298.         //     [ Cdot ] = [ dAcc/dPos ] * [ A ] + [ dAcc/dVel ] * [ C ]
  299.         //     [ Ddot ] = [ dAcc/dPos ] * [ B ] + [ dAcc/dVel ] * [ D ]

  300.         // The following loops compute these expressions taking care of the mapping of the
  301.         // (A, B, C, D) matrices into the single dimension array p and of the mapping of the
  302.         // (Adot, Bdot, Cdot, Ddot) matrices into the single dimension array pDot.

  303.         // copy C and E into Adot and Bdot
  304.         final int stateDim = 6;
  305.         final double[] p = s.getAdditionalState(getName());
  306.         final double[] pDot = new double[p.length];
  307.         System.arraycopy(p, dimEpoch * stateDim, pDot, 0, dimEpoch * stateDim);

  308.         // compute Cdot and Ddot
  309.         for (int i = 0; i < dimEpoch; ++i) {
  310.             final double[] dAdPi = dAccdPos[i];
  311.             final double[] dAdVi = dAccdVel[i];
  312.             for (int j = 0; j < stateDim; ++j) {
  313.                 pDot[(dimEpoch + i) * stateDim + j] =
  314.                     dAdPi[0] * p[j]                + dAdPi[1] * p[j +     stateDim] + dAdPi[2] * p[j + 2 * stateDim] +
  315.                     dAdVi[0] * p[j + 3 * stateDim] + dAdVi[1] * p[j + 4 * stateDim] + dAdVi[2] * p[j + 5 * stateDim];
  316.             }
  317.         }

  318.         for (int k = 0; k < paramDimEpoch; ++k) {
  319.             // the variational equations of the parameters Jacobian matrix are computed
  320.             // one column at a time, they have the following form:
  321.             // [      ]   [                 |                  ]   [   ]   [                  ]
  322.             // [ Edot ]   [  dVel/dPos = 0  |  dVel/dVel = Id  ]   [ E ]   [  dVel/dParam = 0 ]
  323.             // [      ]   [                 |                  ]   [   ]   [                  ]
  324.             // --------   ------------------+------------------- * ----- + --------------------
  325.             // [      ]   [                 |                  ]   [   ]   [                  ]
  326.             // [ Fdot ] = [    dAcc/dPos    |     dAcc/dVel    ]   [ F ]   [    dAcc/dParam   ]
  327.             // [      ]   [                 |                  ]   [   ]   [                  ]

  328.             // The E and F sub-columns and their derivatives (Edot, Fdot) are 3 elements columns.

  329.             // The expanded multiplication and addition above can be rewritten to take into
  330.             // account the fixed values found in the sub-matrices in the left factor. This leads to:

  331.             //     [ Edot ] = [ F ]
  332.             //     [ Fdot ] = [ dAcc/dPos ] * [ E ] + [ dAcc/dVel ] * [ F ] + [ dAcc/dParam ]

  333.             // The following loops compute these expressions taking care of the mapping of the
  334.             // (E, F) columns into the single dimension array p and of the mapping of the
  335.             // (Edot, Fdot) columns into the single dimension array pDot.

  336.             // copy F into Edot
  337.             final int columnTop = stateDim * stateDim + k;
  338.             pDot[columnTop]                     = p[columnTop + 3 * paramDimEpoch];
  339.             pDot[columnTop +     paramDimEpoch] = p[columnTop + 4 * paramDimEpoch];
  340.             pDot[columnTop + 2 * paramDimEpoch] = p[columnTop + 5 * paramDimEpoch];

  341.             // compute Fdot
  342.             for (int i = 0; i < dimEpoch; ++i) {
  343.                 final double[] dAdP = dAccdPos[i];
  344.                 final double[] dAdV = dAccdVel[i];
  345.                 pDot[columnTop + (dimEpoch + i) * paramDimEpoch] =
  346.                     dAccdParam[i][k] +
  347.                     dAdP[0] * p[columnTop]                     + dAdP[1] * p[columnTop +     paramDimEpoch] + dAdP[2] * p[columnTop + 2 * paramDimEpoch] +
  348.                     dAdV[0] * p[columnTop + 3 * paramDimEpoch] + dAdV[1] * p[columnTop + 4 * paramDimEpoch] + dAdV[2] * p[columnTop + 5 * paramDimEpoch];
  349.             }

  350.         }

  351.         return pDot;

  352.     }

  353.     /** Fill Jacobians rows.
  354.      * @param derivatives derivatives of a component of acceleration (along either x, y or z)
  355.      * @param index component index (0 for x, 1 for y, 2 for z)
  356.      * @param freeStateParameters number of free parameters, either 3 (position),
  357.      * 6 (position-velocity) or 7 (position-velocity-mass)
  358.      * @param dAccdPos Jacobian of acceleration with respect to spacecraft position
  359.      * @param dAccdVel Jacobian of acceleration with respect to spacecraft velocity
  360.      */
  361.     private void addToRow(final double[] derivatives, final int index, final int freeStateParameters,
  362.                           final double[][] dAccdPos, final double[][] dAccdVel) {

  363.         for (int i = 0; i < 3; ++i) {
  364.             dAccdPos[index][i] += derivatives[i];
  365.         }
  366.         if (freeStateParameters > 3) {
  367.             for (int i = 0; i < 3; ++i) {
  368.                 dAccdVel[index][i] += derivatives[i + 3];
  369.             }
  370.         }

  371.     }

  372. }