NewtonFixedBoundaryCartesianSingleShooting.java

  1. /* Copyright 2022-2025 Romain Serra
  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.control.indirect.shooting;

  18. import org.hipparchus.analysis.differentiation.Gradient;
  19. import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
  20. import org.hipparchus.geometry.euclidean.threed.Vector3D;
  21. import org.hipparchus.linear.DecompositionSolver;
  22. import org.hipparchus.linear.LUDecomposition;
  23. import org.hipparchus.linear.MatrixUtils;
  24. import org.hipparchus.linear.RealMatrix;
  25. import org.hipparchus.linear.RealVector;
  26. import org.hipparchus.util.FastMath;
  27. import org.orekit.control.indirect.shooting.boundary.CartesianBoundaryConditionChecker;
  28. import org.orekit.control.indirect.shooting.boundary.FixedTimeBoundaryOrbits;
  29. import org.orekit.control.indirect.shooting.boundary.FixedTimeCartesianBoundaryStates;
  30. import org.orekit.control.indirect.shooting.propagation.ShootingPropagationSettings;
  31. import org.orekit.propagation.FieldSpacecraftState;
  32. import org.orekit.utils.FieldPVCoordinates;

  33. /**
  34.  * Class for indirect single shooting methods with Cartesian coordinates for fixed time fixed boundary.
  35.  * Update is the classical Newton-Raphson one. It is computed using an LU matrix decomposition.
  36.  *
  37.  * @author Romain Serra
  38.  * @since 12.2
  39.  */
  40. public class NewtonFixedBoundaryCartesianSingleShooting extends AbstractFixedBoundaryCartesianSingleShooting {

  41.     /** Default value for singularity threshold. */
  42.     private static final double DEFAULT_SINGULARITY_THRESHOLD = 1e-11;

  43.     /** Threshold for singularity exception in linear system inversion. */
  44.     private double singularityThreshold = DEFAULT_SINGULARITY_THRESHOLD;

  45.     /** Multiplying (positive) factor for the Newton correction step. */
  46.     private double stepFactor = 1.;

  47.     /**
  48.      * Constructor with boundary conditions as orbits.
  49.      * @param propagationSettings propagation settings
  50.      * @param boundaryConditions boundary conditions as {@link FixedTimeCartesianBoundaryStates}
  51.      * @param convergenceChecker convergence checker
  52.      */
  53.     public NewtonFixedBoundaryCartesianSingleShooting(final ShootingPropagationSettings propagationSettings,
  54.                                                       final FixedTimeCartesianBoundaryStates boundaryConditions,
  55.                                                       final CartesianBoundaryConditionChecker convergenceChecker) {
  56.         super(propagationSettings, boundaryConditions, convergenceChecker);
  57.     }

  58.     /**
  59.      * Constructor with boundary conditions as orbits.
  60.      * @param propagationSettings propagation settings
  61.      * @param boundaryConditions boundary conditions as {@link FixedTimeBoundaryOrbits}
  62.      * @param convergenceChecker convergence checker
  63.      */
  64.     public NewtonFixedBoundaryCartesianSingleShooting(final ShootingPropagationSettings propagationSettings,
  65.                                                       final FixedTimeBoundaryOrbits boundaryConditions,
  66.                                                       final CartesianBoundaryConditionChecker convergenceChecker) {
  67.         super(propagationSettings, boundaryConditions, convergenceChecker);
  68.     }

  69.     @Override
  70.     public int getMaximumIterationCount() {
  71.         return getConditionChecker().getMaximumIterationCount();
  72.     }

  73.     /**
  74.      * Setter for singularity threshold in LU decomposition.
  75.      * @param singularityThreshold new threshold value
  76.      * @since 13.0
  77.      */
  78.     public void setSingularityThreshold(final double singularityThreshold) {
  79.         this.singularityThreshold = singularityThreshold;
  80.     }

  81.     /**
  82.      * Getter for singularity threshold in LU decomposition.
  83.      * @return threshold
  84.      * @since 13.0
  85.      */
  86.     public double getSingularityThreshold() {
  87.         return singularityThreshold;
  88.     }

  89.     /**
  90.      * Setter for the step factor.
  91.      * @param stepFactor new value for the step factor
  92.      * @since 13.0
  93.      */
  94.     public void setStepFactor(final double stepFactor) {
  95.         this.stepFactor = FastMath.abs(stepFactor);
  96.     }

  97.     /** {@inheritDoc} */
  98.     @Override
  99.     protected double[] updateShootingVariables(final double[] originalShootingVariables,
  100.                                                final FieldSpacecraftState<Gradient> fieldTerminalState) {
  101.         // form defects and their Jacobian matrix
  102.         final double[] defects = new double[originalShootingVariables.length];
  103.         final double[][] defectsJacobianData = new double[defects.length][defects.length];
  104.         final double reciprocalScalePosition = 1. / getScalePositionDefects();
  105.         final double reciprocalScaleVelocity = 1. / getScaleVelocityDefects();
  106.         final FieldPVCoordinates<Gradient> terminalPV = fieldTerminalState.getPVCoordinates();
  107.         final FieldVector3D<Gradient> fieldScaledTerminalPosition = terminalPV.getPosition().scalarMultiply(reciprocalScalePosition);
  108.         final FieldVector3D<Gradient> fieldScaledTerminalVelocity = terminalPV.getVelocity().scalarMultiply(reciprocalScaleVelocity);
  109.         final Vector3D terminalScaledPosition = fieldScaledTerminalPosition.toVector3D();
  110.         final Vector3D terminalScaledVelocity = fieldScaledTerminalVelocity.toVector3D();
  111.         final Vector3D targetScaledPosition = getTerminalCartesianState().getPosition().scalarMultiply(reciprocalScalePosition);
  112.         final Vector3D targetScaledVelocity = getTerminalCartesianState().getVelocity().scalarMultiply(reciprocalScaleVelocity);
  113.         defects[0] = terminalScaledPosition.getX() - targetScaledPosition.getX();
  114.         defectsJacobianData[0] = fieldScaledTerminalPosition.getX().getGradient();
  115.         defects[1] = terminalScaledPosition.getY() - targetScaledPosition.getY();
  116.         defectsJacobianData[1] = fieldScaledTerminalPosition.getY().getGradient();
  117.         defects[2] = terminalScaledPosition.getZ() - targetScaledPosition.getZ();
  118.         defectsJacobianData[2] = fieldScaledTerminalPosition.getZ().getGradient();
  119.         defects[3] = terminalScaledVelocity.getX() - targetScaledVelocity.getX();
  120.         defectsJacobianData[3] = fieldScaledTerminalVelocity.getX().getGradient();
  121.         defects[4] = terminalScaledVelocity.getY() - targetScaledVelocity.getY();
  122.         defectsJacobianData[4] = fieldScaledTerminalVelocity.getY().getGradient();
  123.         defects[5] = terminalScaledVelocity.getZ() - targetScaledVelocity.getZ();
  124.         defectsJacobianData[5] = fieldScaledTerminalVelocity.getZ().getGradient();
  125.         if (originalShootingVariables.length != 6) {
  126.             final String adjointName = getPropagationSettings().getAdjointDynamicsProvider().getAdjointName();
  127.             final Gradient terminalMassAdjoint = fieldTerminalState.getAdditionalState(adjointName)[6];
  128.             defects[6] = terminalMassAdjoint.getValue();
  129.             defectsJacobianData[6] = terminalMassAdjoint.getGradient();
  130.         }
  131.         // apply Newton's formula
  132.         final double[] correction = computeCorrection(defects, defectsJacobianData);
  133.         final double[] correctedAdjoint = originalShootingVariables.clone();
  134.         for (int i = 0; i < correction.length; i++) {
  135.             correctedAdjoint[i] += correction[i];
  136.         }
  137.         return correctedAdjoint;
  138.     }

  139.     /**
  140.      * Compute Newton-Raphson correction.
  141.      * @param defects defects
  142.      * @param defectsJacobianData Jacobian matrix of defects w.r.t. shooting variables
  143.      * @return correction to shooting variables
  144.      */
  145.     private double[] computeCorrection(final double[] defects, final double[][] defectsJacobianData) {
  146.         final RealMatrix defectsJacobian = MatrixUtils.createRealMatrix(defectsJacobianData);
  147.         final DecompositionSolver solver = new LUDecomposition(defectsJacobian, singularityThreshold).getSolver();
  148.         final RealVector negatedDefects = MatrixUtils.createRealVector(defects).mapMultiply(-1);
  149.         final RealVector solved = solver.solve(negatedDefects);
  150.         final double[] corrections = new double[solved.getDimension()];
  151.         for (int i = 0; i < corrections.length; i++) {
  152.             corrections[i] = solved.getEntry(i) * getScales()[i] * stepFactor;
  153.         }
  154.         return corrections;
  155.     }
  156. }