AbstractAnalyticalMatricesHarvester.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.analytical;

  18. import java.util.Arrays;
  19. import java.util.List;

  20. import org.hipparchus.analysis.differentiation.Gradient;
  21. import org.hipparchus.linear.MatrixUtils;
  22. import org.hipparchus.linear.RealMatrix;
  23. import org.orekit.orbits.FieldOrbit;
  24. import org.orekit.propagation.AbstractMatricesHarvester;
  25. import org.orekit.propagation.AdditionalStateProvider;
  26. import org.orekit.propagation.FieldSpacecraftState;
  27. import org.orekit.propagation.SpacecraftState;
  28. import org.orekit.time.AbsoluteDate;
  29. import org.orekit.time.FieldAbsoluteDate;
  30. import org.orekit.utils.DoubleArrayDictionary;
  31. import org.orekit.utils.FieldPVCoordinates;
  32. import org.orekit.utils.ParameterDriver;

  33. /**
  34.  * Base class harvester between two-dimensional Jacobian
  35.  * matrices and analytical orbit propagator.
  36.  * @author Thomas Paulet
  37.  * @author Bryan Cazabonne
  38.  * @since 11.1
  39.  */
  40. public abstract class AbstractAnalyticalMatricesHarvester extends AbstractMatricesHarvester implements AdditionalStateProvider {

  41.     /** Columns names for parameters. */
  42.     private List<String> columnsNames;

  43.     /** Epoch of the last computed state transition matrix. */
  44.     private AbsoluteDate epoch;

  45.     /** Analytical derivatives that apply to State Transition Matrix. */
  46.     private final double[][] analyticalDerivativesStm;

  47.     /** Analytical derivatives that apply to Jacobians columns. */
  48.     private final DoubleArrayDictionary analyticalDerivativesJacobianColumns;

  49.     /** Propagator bound to this harvester. */
  50.     private final AbstractAnalyticalPropagator propagator;

  51.     /** Simple constructor.
  52.      * <p>
  53.      * The arguments for initial matrices <em>must</em> be compatible with the
  54.      * {@link org.orekit.orbits.OrbitType orbit type}
  55.      * and {@link org.orekit.orbits.PositionAngle position angle} that will be used by propagator
  56.      * </p>
  57.      * @param propagator propagator bound to this harvester
  58.      * @param stmName State Transition Matrix state name
  59.      * @param initialStm initial State Transition Matrix ∂Y/∂Y₀,
  60.      * if null (which is the most frequent case), assumed to be 6x6 identity
  61.      * @param initialJacobianColumns initial columns of the Jacobians matrix with respect to parameters,
  62.      * if null or if some selected parameters are missing from the dictionary, the corresponding
  63.      * initial column is assumed to be 0
  64.      */
  65.     protected AbstractAnalyticalMatricesHarvester(final AbstractAnalyticalPropagator propagator, final String stmName,
  66.                                                   final RealMatrix initialStm, final DoubleArrayDictionary initialJacobianColumns) {
  67.         super(stmName, initialStm, initialJacobianColumns);
  68.         this.propagator                           = propagator;
  69.         this.epoch                                = propagator.getInitialState().getDate();
  70.         this.columnsNames                         = null;
  71.         this.analyticalDerivativesStm             = getInitialStateTransitionMatrix().getData();
  72.         this.analyticalDerivativesJacobianColumns = new DoubleArrayDictionary();
  73.     }

  74.     /** {@inheritDoc} */
  75.     @Override
  76.     public List<String> getJacobiansColumnsNames() {
  77.         return columnsNames == null ? propagator.getJacobiansColumnsNames() : columnsNames;
  78.     }

  79.     /** {@inheritDoc} */
  80.     @Override
  81.     public void freezeColumnsNames() {
  82.         columnsNames = getJacobiansColumnsNames();
  83.     }

  84.     /** {@inheritDoc} */
  85.     @Override
  86.     public String getName() {
  87.         return getStmName();
  88.     }

  89.     /** {@inheritDoc} */
  90.     @Override
  91.     public double[] getAdditionalState(final SpacecraftState state) {
  92.         // Update the partial derivatives if needed
  93.         updateDerivativesIfNeeded(state);
  94.         // Return the state transition matrix in an array
  95.         return toArray(analyticalDerivativesStm);
  96.     }

  97.     /** {@inheritDoc} */
  98.     @Override
  99.     public RealMatrix getStateTransitionMatrix(final SpacecraftState state) {
  100.         // Check if additional state is defined
  101.         if (!state.hasAdditionalState(getName())) {
  102.             return null;
  103.         }
  104.         // Return the state transition matrix
  105.         return toRealMatrix(state.getAdditionalState(getName()));
  106.     }

  107.     /** {@inheritDoc} */
  108.     @Override
  109.     public RealMatrix getParametersJacobian(final SpacecraftState state) {
  110.         // Update the partial derivatives if needed
  111.         updateDerivativesIfNeeded(state);

  112.         // Estimated parameters
  113.         final List<String> names = getJacobiansColumnsNames();
  114.         if (names == null || names.isEmpty()) {
  115.             return null;
  116.         }

  117.         // Initialize Jacobian
  118.         final RealMatrix dYdP = MatrixUtils.createRealMatrix(STATE_DIMENSION, names.size());

  119.         // Add the derivatives
  120.         for (int j = 0; j < names.size(); ++j) {
  121.             final double[] column = analyticalDerivativesJacobianColumns.get(names.get(j));
  122.             if (column != null) {
  123.                 for (int i = 0; i < STATE_DIMENSION; i++) {
  124.                     dYdP.addToEntry(i, j, column[i]);
  125.                 }
  126.             }
  127.         }

  128.         // Return
  129.         return dYdP;
  130.     }

  131.     /** {@inheritDoc} */
  132.     @Override
  133.     public void setReferenceState(final SpacecraftState reference) {

  134.         // reset derivatives to zero
  135.         for (final double[] row : analyticalDerivativesStm) {
  136.             Arrays.fill(row, 0.0);
  137.         }
  138.         analyticalDerivativesJacobianColumns.clear();

  139.         final AbstractAnalyticalGradientConverter converter           = getGradientConverter();
  140.         final FieldSpacecraftState<Gradient> gState                   = converter.getState();
  141.         final Gradient[] gParameters                                  = converter.getParameters(gState);
  142.         final FieldAbstractAnalyticalPropagator<Gradient> gPropagator = converter.getPropagator(gState, gParameters);

  143.         // Compute Jacobian
  144.         final AbsoluteDate target               = reference.getDate();
  145.         final FieldAbsoluteDate<Gradient> start = gPropagator.getInitialState().getDate();
  146.         final double dt                         = target.durationFrom(start.toAbsoluteDate());
  147.         final FieldOrbit<Gradient> gOrbit       = gPropagator.propagateOrbit(start.shiftedBy(dt), gParameters);
  148.         final FieldPVCoordinates<Gradient> gPv  = gOrbit.getPVCoordinates();

  149.         final double[] derivativesX   = gPv.getPosition().getX().getGradient();
  150.         final double[] derivativesY   = gPv.getPosition().getY().getGradient();
  151.         final double[] derivativesZ   = gPv.getPosition().getZ().getGradient();
  152.         final double[] derivativesVx  = gPv.getVelocity().getX().getGradient();
  153.         final double[] derivativesVy  = gPv.getVelocity().getY().getGradient();
  154.         final double[] derivativesVz  = gPv.getVelocity().getZ().getGradient();

  155.         // Update Jacobian with respect to state
  156.         addToRow(derivativesX,  0);
  157.         addToRow(derivativesY,  1);
  158.         addToRow(derivativesZ,  2);
  159.         addToRow(derivativesVx, 3);
  160.         addToRow(derivativesVy, 4);
  161.         addToRow(derivativesVz, 5);

  162.         // Partial derivatives of the state with respect to propagation parameters
  163.         int paramsIndex = converter.getFreeStateParameters();
  164.         for (ParameterDriver driver : converter.getParametersDrivers()) {
  165.             if (driver.isSelected()) {

  166.                 // get the partials derivatives for this driver
  167.                 DoubleArrayDictionary.Entry entry = analyticalDerivativesJacobianColumns.getEntry(driver.getName());
  168.                 if (entry == null) {
  169.                     // create an entry filled with zeroes
  170.                     analyticalDerivativesJacobianColumns.put(driver.getName(), new double[STATE_DIMENSION]);
  171.                     entry = analyticalDerivativesJacobianColumns.getEntry(driver.getName());
  172.                 }

  173.                 // add the contribution of the current force model
  174.                 entry.increment(new double[] {
  175.                     derivativesX[paramsIndex], derivativesY[paramsIndex], derivativesZ[paramsIndex],
  176.                     derivativesVx[paramsIndex], derivativesVy[paramsIndex], derivativesVz[paramsIndex]
  177.                 });
  178.                 ++paramsIndex;

  179.             }
  180.         }

  181.         // Update the epoch of the last computed partial derivatives
  182.         epoch = target;

  183.     }

  184.     /** Update the partial derivatives (if needed).
  185.      * @param state current spacecraft state
  186.      */
  187.     private void updateDerivativesIfNeeded(final SpacecraftState state) {
  188.         if (!state.getDate().isEqualTo(epoch)) {
  189.             setReferenceState(state);
  190.         }
  191.     }

  192.     /** Fill State Transition Matrix rows.
  193.      * @param derivatives derivatives of a component
  194.      * @param index component index
  195.      */
  196.     private void addToRow(final double[] derivatives, final int index) {
  197.         for (int i = 0; i < 6; i++) {
  198.             analyticalDerivativesStm[index][i] += derivatives[i];
  199.         }
  200.     }

  201.     /** Convert an array to a matrix (6x6 dimension).
  202.      * @param array input array
  203.      * @return the corresponding matrix
  204.      */
  205.     private RealMatrix toRealMatrix(final double[] array) {
  206.         final RealMatrix matrix = MatrixUtils.createRealMatrix(STATE_DIMENSION, STATE_DIMENSION);
  207.         int index = 0;
  208.         for (int i = 0; i < STATE_DIMENSION; ++i) {
  209.             for (int j = 0; j < STATE_DIMENSION; ++j) {
  210.                 matrix.setEntry(i, j, array[index++]);
  211.             }
  212.         }
  213.         return matrix;
  214.     }

  215.     /** Set the STM data into an array.
  216.      * @param matrix STM matrix
  217.      * @return an array containing the STM data
  218.      */
  219.     private double[] toArray(final double[][] matrix) {
  220.         final double[] array = new double[STATE_DIMENSION * STATE_DIMENSION];
  221.         int index = 0;
  222.         for (int i = 0; i < STATE_DIMENSION; ++i) {
  223.             final double[] row = matrix[i];
  224.             for (int j = 0; j < STATE_DIMENSION; ++j) {
  225.                 array[index++] = row[j];
  226.             }
  227.         }
  228.         return array;
  229.     }

  230.     /**
  231.      * Get the gradient converter related to the analytical orbit propagator.
  232.      * @return the gradient converter
  233.      */
  234.     public abstract AbstractAnalyticalGradientConverter getGradientConverter();

  235. }