IntegrableJacobianColumnGenerator.java

  1. /* Copyright 2002-2024 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 org.orekit.propagation.SpacecraftState;
  19. import org.orekit.propagation.integration.AdditionalDerivativesProvider;
  20. import org.orekit.propagation.integration.CombinedDerivatives;

  21. /** Generator for one column of a Jacobian matrix.
  22.  * <p>
  23.  * This generator is based on variational equations, so
  24.  * it implements {@link AdditionalDerivativesProvider} and
  25.  * computes only the derivative of the Jacobian column, to
  26.  * be integrated by the propagator alongside the primary state.
  27.  * </p>
  28.  * @author Luc Maisonobe
  29.  * @since 11.1
  30.  */
  31. class IntegrableJacobianColumnGenerator
  32.     implements AdditionalDerivativesProvider, StateTransitionMatrixGenerator.PartialsObserver {

  33.     /** Name of the state for State Transition Matrix. */
  34.     private final String stmName;

  35.     /** Name of the parameter corresponding to the column. */
  36.     private final String columnName;

  37.     /** Last value computed for the partial derivatives. */
  38.     private final double[] pDot;

  39.     /** Simple constructor.
  40.      * <p>
  41.      * The generator for State Transition Matrix <em>must</em> be registered as
  42.      * an integrable generator to the same propagator as the instance, as it
  43.      * must be scheduled to update the state before the instance
  44.      * </p>
  45.      * @param stmGenerator generator for State Transition Matrix
  46.      * @param columnName name of the parameter corresponding to the column
  47.      */
  48.     IntegrableJacobianColumnGenerator(final StateTransitionMatrixGenerator stmGenerator, final String columnName) {
  49.         this.stmName    = stmGenerator.getName();
  50.         this.columnName = columnName;
  51.         this.pDot       = new double[getDimension()];
  52.         stmGenerator.addObserver(columnName, this);
  53.     }

  54.     /** {@inheritDoc} */
  55.     @Override
  56.     public String getName() {
  57.         return columnName;
  58.     }

  59.     /** Get the dimension of the generated column.
  60.      * @return dimension of the generated column
  61.      */
  62.     public int getDimension() {
  63.         return 6;
  64.     }

  65.     /** {@inheritDoc}
  66.      * <p>
  67.      * The column derivative can be computed only if the State Transition Matrix derivatives
  68.      * are available, as it implies the STM generator has already been run.
  69.      * </p>
  70.      */
  71.     @Override
  72.     public boolean yields(final SpacecraftState state) {
  73.         return !state.hasAdditionalStateDerivative(stmName);
  74.     }

  75.     /** {@inheritDoc} */
  76.     @Override
  77.     public void partialsComputed(final SpacecraftState state, final double[] factor, final double[] accelerationPartials) {
  78.         // retrieve current Jacobian column
  79.         final double[] p = state.getAdditionalState(getName());

  80.         // compute time derivative of the Jacobian column
  81.         StateTransitionMatrixGenerator.multiplyMatrix(factor, p, pDot, 1);
  82.         pDot[3] += accelerationPartials[0];
  83.         pDot[4] += accelerationPartials[1];
  84.         pDot[5] += accelerationPartials[2];
  85.     }

  86.     /** {@inheritDoc} */
  87.     @Override
  88.     public CombinedDerivatives combinedDerivatives(final SpacecraftState s) {
  89.         return new CombinedDerivatives(pDot, null);
  90.     }

  91. }