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.propagation.numerical;
18  
19  import org.hipparchus.geometry.euclidean.threed.Vector3D;
20  import org.hipparchus.util.FastMath;
21  import org.orekit.errors.OrekitIllegalArgumentException;
22  import org.orekit.errors.OrekitMessages;
23  import org.orekit.forces.ForceModel;
24  import org.orekit.orbits.OrbitType;
25  import org.orekit.orbits.PositionAngleType;
26  import org.orekit.propagation.SpacecraftState;
27  
28  import java.util.Arrays;
29  import java.util.List;
30  
31  /**
32   * First-order differential model representing the equations of motion.
33   *
34   * @author Luc Maisonobe
35   * @author Romain Serra
36   * @since 13.1
37   * @see NumericalPropagator
38   */
39  class NumericalTimeDerivativesEquations implements TimeDerivativesEquations {
40  
41      /** Derivatives array. */
42      private final double[] yDot;
43  
44      /** Orbit type of coordinates. */
45      private final OrbitType orbitType;
46  
47      /** Position angle type of coordinates. */
48      private final PositionAngleType positionAngleType;
49  
50      /** Forces list. */
51      private final List<ForceModel> forceModels;
52  
53      /** Jacobian of the orbital parameters with respect to the Cartesian parameters. */
54      private double[][] coordinatesJacobian;
55  
56      /** Current state. */
57      private SpacecraftState currentState;
58  
59      /**
60       * Constructor.
61       * @param orbitType orbit type used in equations
62       * @param positionAngleType angle type used in equations if applicable
63       * @param forceModels forces
64       */
65      NumericalTimeDerivativesEquations(final OrbitType orbitType, final PositionAngleType positionAngleType,
66                                        final List<ForceModel> forceModels) {
67          this.orbitType = orbitType;
68          this.positionAngleType = positionAngleType;
69          this.forceModels = forceModels;
70          this.yDot     = new double[7];
71          this.coordinatesJacobian = new double[6][6];
72          // default value for Jacobian is identity
73          for (int i = 0; i < coordinatesJacobian.length; ++i) {
74              Arrays.fill(coordinatesJacobian[i], 0.0);
75              coordinatesJacobian[i][i] = 1.0;
76          }
77      }
78  
79      /**
80       * Getter for the force models.
81       * @return forces
82       */
83      List<ForceModel> getForceModels() {
84          return forceModels;
85      }
86  
87      /**
88       * Setter for the coordinates' Jacobian matrix.
89       * @param coordinatesJacobian matrix
90       */
91      void setCoordinatesJacobian(final double[][] coordinatesJacobian) {
92          this.coordinatesJacobian = coordinatesJacobian;
93      }
94  
95      /**
96       * Setter for the current state.
97       * @param currentState state
98       */
99      void setCurrentState(final SpacecraftState currentState) {
100         this.currentState = currentState;
101     }
102 
103     /**
104      * Compute first-order, time derivatives a.k.a. rates.
105      * @param state state where to evaluate derivatives
106      * @return derivatives
107      */
108     double[] computeTimeDerivatives(final SpacecraftState state) {
109         currentState = state;
110         Arrays.fill(yDot, 0.0);
111 
112         // compute the contributions of all perturbing forces,
113         // using the Kepler contribution at the end since
114         // NewtonianAttraction is always the last instance in the list
115         for (final ForceModel forceModel : forceModels) {
116             forceModel.addContribution(state, this);
117         }
118 
119         if (orbitType == null) {
120             // position derivative is velocity, and was not added above in the force models
121             // (it is added when orbit type is non-null because NewtonianAttraction considers it)
122             final Vector3D velocity = state.getPVCoordinates().getVelocity();
123             yDot[0] += velocity.getX();
124             yDot[1] += velocity.getY();
125             yDot[2] += velocity.getZ();
126         }
127 
128         return yDot.clone();
129     }
130 
131     /** {@inheritDoc} */
132     @Override
133     public void addKeplerContribution(final double mu) {
134         if (orbitType == null) {
135             // if mu is neither 0 nor NaN, we want to include Newtonian acceleration
136             if (mu > 0) {
137                 // velocity derivative is Newtonian acceleration
138                 final Vector3D position = currentState.getPosition();
139                 final double r2         = position.getNormSq();
140                 final double coeff      = -mu / (r2 * FastMath.sqrt(r2));
141                 yDot[3] += coeff * position.getX();
142                 yDot[4] += coeff * position.getY();
143                 yDot[5] += coeff * position.getZ();
144             }
145 
146         } else {
147             // propagation uses regular orbits
148             orbitType.convertType(currentState.getOrbit()).addKeplerContribution(positionAngleType, mu, yDot);
149         }
150     }
151 
152     /** {@inheritDoc} */
153     @Override
154     public void addNonKeplerianAcceleration(final Vector3D gamma) {
155         for (int i = 0; i < 6; ++i) {
156             final double[] jRow = coordinatesJacobian[i];
157             yDot[i] += jRow[3] * gamma.getX() + jRow[4] * gamma.getY() + jRow[5] * gamma.getZ();
158         }
159     }
160 
161     /** {@inheritDoc} */
162     @Override
163     public void addMassDerivative(final double q) {
164         if (q > 0) {
165             throw new OrekitIllegalArgumentException(OrekitMessages.POSITIVE_FLOW_RATE, q);
166         }
167         yDot[6] += q;
168     }
169 }