STMEquations.java

  1. /* Copyright 2002-2020 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.cr3bp;

  18. import java.util.Arrays;

  19. import org.hipparchus.analysis.differentiation.DerivativeStructure;
  20. import org.hipparchus.linear.Array2DRowRealMatrix;
  21. import org.hipparchus.linear.RealMatrix;
  22. import org.orekit.bodies.CR3BPSystem;
  23. import org.orekit.propagation.SpacecraftState;
  24. import org.orekit.propagation.integration.AdditionalEquations;

  25. /** Class calculating the state transition matrix coefficient for CR3BP Computation.
  26.  * @see "Dynamical systems, the three-body problem, and space mission design, Koon, Lo, Marsden, Ross"
  27.  * @author Vincent Mouraux
  28.  * @since 10.2
  29.  */
  30. public class STMEquations implements AdditionalEquations {

  31.     /** Matrix Dimension. */
  32.     private static final int DIM = 6;

  33.     /** Mass ratio of the considered CR3BP System. */
  34.     private final CR3BPSystem syst;

  35.     /** Name of the equations. */
  36.     private final String name;

  37.     /** Potential Hessian Matrix. */
  38.     private final double[][] jacobian = new double[DIM][DIM];

  39.     /** Simple constructor.
  40.      * @param syst CR3BP System considered
  41.      */
  42.     public STMEquations(final CR3BPSystem syst) {
  43.         this.syst = syst;
  44.         this.name = "stmEquations";

  45.         // Jacobian constant values initialization
  46.         for (int j = 0; j < jacobian.length; ++j) {
  47.             Arrays.fill(jacobian[j], 0.0);
  48.         }

  49.         jacobian[0][3] = 1.0;
  50.         jacobian[1][4] = 1.0;
  51.         jacobian[2][5] = 1.0;
  52.         jacobian[3][4] = 2.0;
  53.         jacobian[4][3] = -2.0;
  54.     }

  55.     /** Method adding the standard initial values of the additional state to the initial spacecraft state.
  56.      * @param s Initial state of the system
  57.      * @return s Initial augmented (with the additional equations) state
  58.      */
  59.     public SpacecraftState setInitialPhi(final SpacecraftState s) {
  60.         final int stateDimension = 36;
  61.         final double[] phi = new double[stateDimension];
  62.         for (int i = 0; i < stateDimension; ++i) {
  63.             phi[i] = 0.0;
  64.         }

  65.         for (int i = 0; i < stateDimension; i = i + 7) {
  66.             phi[i] = 1.0;
  67.         }

  68.         return s.addAdditionalState(name, phi);
  69.     }

  70.     /** {@inheritDoc} */
  71.     public double[] computeDerivatives(final SpacecraftState s, final double[] dPhi) {

  72.         // State Transition Matrix
  73.         final double[] phi = s.getAdditionalState(getName());

  74.         // Spacecraft Potential
  75.         final DerivativeStructure potential = new CR3BPForceModel(syst).getPotential(s);

  76.         // Potential derivatives
  77.         final double[] dU = potential.getAllDerivatives();

  78.         // second order derivatives index
  79.         final int idXX = potential.getFactory().getCompiler().getPartialDerivativeIndex(2, 0, 0);
  80.         final int idXY = potential.getFactory().getCompiler().getPartialDerivativeIndex(1, 1, 0);
  81.         final int idXZ = potential.getFactory().getCompiler().getPartialDerivativeIndex(1, 0, 1);
  82.         final int idYY = potential.getFactory().getCompiler().getPartialDerivativeIndex(0, 2, 0);
  83.         final int idYZ = potential.getFactory().getCompiler().getPartialDerivativeIndex(0, 1, 1);
  84.         final int idZZ = potential.getFactory().getCompiler().getPartialDerivativeIndex(0, 0, 2);

  85.         // New Jacobian values
  86.         jacobian[3][0] = dU[idXX];
  87.         jacobian[4][1] = dU[idYY];
  88.         jacobian[5][2] = dU[idZZ];
  89.         jacobian[3][1] = dU[idXY];
  90.         jacobian[4][0] = jacobian[3][1];
  91.         jacobian[3][2] = dU[idXZ];
  92.         jacobian[5][0] = jacobian[3][2];
  93.         jacobian[4][2] = dU[idYZ];
  94.         jacobian[5][1] = jacobian[4][2];

  95.         // STM derivatives computation : dPhi = Jacobian * Phi if both dPhi and Phi are defined as Matrix
  96.         for (int k = 0; k < DIM; k++) {
  97.             for (int l = 0; l < DIM; l++) {
  98.                 for (int i = 0; i < DIM; i++) {
  99.                     dPhi[DIM * k + l] =
  100.                         dPhi[DIM * k + l] + jacobian[k][i] * phi[DIM * i + l];
  101.                 }
  102.             }
  103.         }

  104.         // these equations have no effect on the main state itself
  105.         return null;
  106.     }

  107.     /** {@inheritDoc} */
  108.     public String getName() {
  109.         return name;
  110.     }

  111.     /** Method returning the State Transition Matrix.
  112.      * @param s SpacecraftState of the system
  113.      * @return phiM State Transition Matrix
  114.      */
  115.     public RealMatrix getStateTransitionMatrix(final SpacecraftState s) {
  116.         final double[][] phi2dA = new double[DIM][DIM];
  117.         final double[] stm = s.getAdditionalState(getName());
  118.         for (int i = 0; i < DIM; i++) {
  119.             for (int j = 0; j < 6; j++) {
  120.                 phi2dA[i][j] = stm[DIM * i + j];
  121.             }
  122.         }
  123.         return new Array2DRowRealMatrix(phi2dA, false);
  124.     }
  125. }