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.FieldSpacecraftState;
  26. import org.orekit.propagation.SpacecraftState;
  27. import org.orekit.time.AbsoluteDate;
  28. import org.orekit.time.FieldAbsoluteDate;
  29. import org.orekit.utils.DoubleArrayDictionary;
  30. import org.orekit.utils.FieldPVCoordinates;
  31. import org.orekit.utils.ParameterDriver;

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

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

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

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

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

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

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

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

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

  83.     /** {@inheritDoc} */
  84.     @Override
  85.     public RealMatrix getStateTransitionMatrix(final SpacecraftState state) {
  86.         // Update the partial derivatives if needed
  87.         updateDerivativesIfNeeded(state);
  88.         //return the state transition matrix
  89.         return MatrixUtils.createRealMatrix(analyticalDerivativesStm);
  90.     }

  91.     /** {@inheritDoc} */
  92.     @Override
  93.     public RealMatrix getParametersJacobian(final SpacecraftState state) {
  94.         // Update the partial derivatives if needed
  95.         updateDerivativesIfNeeded(state);

  96.         // Estimated parameters
  97.         final List<String> names = getJacobiansColumnsNames();
  98.         if (names == null || names.isEmpty()) {
  99.             return null;
  100.         }

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

  103.         // Add the derivatives
  104.         for (int j = 0; j < names.size(); ++j) {
  105.             final double[] column = analyticalDerivativesJacobianColumns.get(names.get(j));
  106.             if (column != null) {
  107.                 for (int i = 0; i < STATE_DIMENSION; i++) {
  108.                     dYdP.addToEntry(i, j, column[i]);
  109.                 }
  110.             }
  111.         }

  112.         // Return
  113.         return dYdP;
  114.     }

  115.     /** {@inheritDoc} */
  116.     @Override
  117.     public void setReferenceState(final SpacecraftState reference) {

  118.         // reset derivatives to zero
  119.         for (final double[] row : analyticalDerivativesStm) {
  120.             Arrays.fill(row, 0.0);
  121.         }
  122.         analyticalDerivativesJacobianColumns.clear();

  123.         final AbstractAnalyticalGradientConverter converter           = getGradientConverter();
  124.         final FieldSpacecraftState<Gradient> gState                   = converter.getState();
  125.         final Gradient[] gParameters                                  = converter.getParameters(gState);
  126.         final FieldAbstractAnalyticalPropagator<Gradient> gPropagator = converter.getPropagator(gState, gParameters);

  127.         // Compute Jacobian
  128.         final AbsoluteDate target               = reference.getDate();
  129.         final FieldAbsoluteDate<Gradient> start = gPropagator.getInitialState().getDate();
  130.         final double dt                         = target.durationFrom(start.toAbsoluteDate());
  131.         final FieldOrbit<Gradient> gOrbit       = gPropagator.propagateOrbit(start.shiftedBy(dt), gParameters);
  132.         final FieldPVCoordinates<Gradient> gPv  = gOrbit.getPVCoordinates();

  133.         final double[] derivativesX   = gPv.getPosition().getX().getGradient();
  134.         final double[] derivativesY   = gPv.getPosition().getY().getGradient();
  135.         final double[] derivativesZ   = gPv.getPosition().getZ().getGradient();
  136.         final double[] derivativesVx  = gPv.getVelocity().getX().getGradient();
  137.         final double[] derivativesVy  = gPv.getVelocity().getY().getGradient();
  138.         final double[] derivativesVz  = gPv.getVelocity().getZ().getGradient();

  139.         // Update Jacobian with respect to state
  140.         addToRow(derivativesX,  0);
  141.         addToRow(derivativesY,  1);
  142.         addToRow(derivativesZ,  2);
  143.         addToRow(derivativesVx, 3);
  144.         addToRow(derivativesVy, 4);
  145.         addToRow(derivativesVz, 5);

  146.         // Partial derivatives of the state with respect to propagation parameters
  147.         int paramsIndex = converter.getFreeStateParameters();
  148.         for (ParameterDriver driver : converter.getParametersDrivers()) {
  149.             if (driver.isSelected()) {

  150.                 // get the partials derivatives for this driver
  151.                 DoubleArrayDictionary.Entry entry = analyticalDerivativesJacobianColumns.getEntry(driver.getName());
  152.                 if (entry == null) {
  153.                     // create an entry filled with zeroes
  154.                     analyticalDerivativesJacobianColumns.put(driver.getName(), new double[STATE_DIMENSION]);
  155.                     entry = analyticalDerivativesJacobianColumns.getEntry(driver.getName());
  156.                 }

  157.                 // add the contribution of the current force model
  158.                 entry.increment(new double[] {
  159.                     derivativesX[paramsIndex], derivativesY[paramsIndex], derivativesZ[paramsIndex],
  160.                     derivativesVx[paramsIndex], derivativesVy[paramsIndex], derivativesVz[paramsIndex]
  161.                 });
  162.                 ++paramsIndex;

  163.             }
  164.         }

  165.         // Update the epoch of the last computed partial derivatives
  166.         epoch = target;

  167.     }

  168.     /** Update the partial derivatives (if needed).
  169.      * @param state current spacecraft state
  170.      */
  171.     private void updateDerivativesIfNeeded(final SpacecraftState state) {
  172.         if (!state.getDate().isEqualTo(epoch)) {
  173.             setReferenceState(state);
  174.         }
  175.     }

  176.     /** Fill State Transition Matrix rows.
  177.      * @param derivatives derivatives of a component
  178.      * @param index component index
  179.      */
  180.     private void addToRow(final double[] derivatives, final int index) {
  181.         for (int i = 0; i < 6; i++) {
  182.             analyticalDerivativesStm[index][i] += derivatives[i];
  183.         }
  184.     }

  185.     /**
  186.      * Get the gradient converter related to the analytical orbit propagator.
  187.      * @return the gradient converter
  188.      */
  189.     public abstract AbstractAnalyticalGradientConverter getGradientConverter();

  190. }