1   /* Copyright 2002-2021 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  
19  import java.util.ArrayList;
20  import java.util.Arrays;
21  import java.util.Collections;
22  import java.util.List;
23  
24  import org.hipparchus.geometry.euclidean.threed.Vector3D;
25  import org.hipparchus.ode.ODEIntegrator;
26  import org.hipparchus.util.FastMath;
27  import org.orekit.annotation.DefaultDataContext;
28  import org.orekit.attitudes.Attitude;
29  import org.orekit.attitudes.AttitudeProvider;
30  import org.orekit.data.DataContext;
31  import org.orekit.errors.OrekitException;
32  import org.orekit.errors.OrekitIllegalArgumentException;
33  import org.orekit.errors.OrekitInternalError;
34  import org.orekit.errors.OrekitMessages;
35  import org.orekit.forces.ForceModel;
36  import org.orekit.forces.gravity.NewtonianAttraction;
37  import org.orekit.forces.inertia.InertialForces;
38  import org.orekit.frames.Frame;
39  import org.orekit.orbits.Orbit;
40  import org.orekit.orbits.OrbitType;
41  import org.orekit.orbits.PositionAngle;
42  import org.orekit.propagation.PropagationType;
43  import org.orekit.propagation.Propagator;
44  import org.orekit.propagation.SpacecraftState;
45  import org.orekit.propagation.events.EventDetector;
46  import org.orekit.propagation.integration.AbstractIntegratedPropagator;
47  import org.orekit.propagation.integration.StateMapper;
48  import org.orekit.time.AbsoluteDate;
49  import org.orekit.utils.AbsolutePVCoordinates;
50  import org.orekit.utils.PVCoordinates;
51  import org.orekit.utils.ParameterDriver;
52  import org.orekit.utils.ParameterObserver;
53  import org.orekit.utils.TimeStampedPVCoordinates;
54  
55  /** This class propagates {@link org.orekit.orbits.Orbit orbits} using
56   * numerical integration.
57   * <p>Numerical propagation is much more accurate than analytical propagation
58   * like for example {@link org.orekit.propagation.analytical.KeplerianPropagator
59   * Keplerian} or {@link org.orekit.propagation.analytical.EcksteinHechlerPropagator
60   * Eckstein-Hechler}, but requires a few more steps to set up to be used properly.
61   * Whereas analytical propagators are configured only thanks to their various
62   * constructors and can be used immediately after construction, numerical propagators
63   * configuration involve setting several parameters between construction time
64   * and propagation time.</p>
65   * <p>The configuration parameters that can be set are:</p>
66   * <ul>
67   *   <li>the initial spacecraft state ({@link #setInitialState(SpacecraftState)})</li>
68   *   <li>the central attraction coefficient ({@link #setMu(double)})</li>
69   *   <li>the various force models ({@link #addForceModel(ForceModel)},
70   *   {@link #removeForceModels()})</li>
71   *   <li>the {@link OrbitType type} of orbital parameters to be used for propagation
72   *   ({@link #setOrbitType(OrbitType)}),
73   *   <li>the {@link PositionAngle type} of position angle to be used in orbital parameters
74   *   to be used for propagation where it is relevant ({@link
75   *   #setPositionAngleType(PositionAngle)}),
76   *   <li>whether {@link org.orekit.propagation.integration.AdditionalEquations additional equations}
77   *   (for example {@link PartialDerivativesEquations Jacobians}) should be propagated along with orbital state
78   *   ({@link #addAdditionalEquations(org.orekit.propagation.integration.AdditionalEquations)}),
79   *   <li>the discrete events that should be triggered during propagation
80   *   ({@link #addEventDetector(EventDetector)},
81   *   {@link #clearEventsDetectors()})</li>
82   *   <li>the binding logic with the rest of the application ({@link #getMultiplexer()})</li>
83   * </ul>
84   * <p>From these configuration parameters, only the initial state is mandatory. The default
85   * propagation settings are in {@link OrbitType#EQUINOCTIAL equinoctial} parameters with
86   * {@link PositionAngle#TRUE true} longitude argument. If the central attraction coefficient
87   * is not explicitly specified, the one used to define the initial orbit will be used.
88   * However, specifying only the initial state and perhaps the central attraction coefficient
89   * would mean the propagator would use only Keplerian forces. In this case, the simpler {@link
90   * org.orekit.propagation.analytical.KeplerianPropagator KeplerianPropagator} class would
91   * perhaps be more effective.</p>
92   * <p>The underlying numerical integrator set up in the constructor may also have its own
93   * configuration parameters. Typical configuration parameters for adaptive stepsize integrators
94   * are the min, max and perhaps start step size as well as the absolute and/or relative errors
95   * thresholds.</p>
96   * <p>The state that is seen by the integrator is a simple seven elements double array.
97   * The six first elements are either:
98   * <ul>
99   *   <li>the {@link org.orekit.orbits.EquinoctialOrbit equinoctial orbit parameters} (a, e<sub>x</sub>,
100  *   e<sub>y</sub>, h<sub>x</sub>, h<sub>y</sub>, λ<sub>M</sub> or λ<sub>E</sub>
101  *   or λ<sub>v</sub>) in meters and radians,</li>
102  *   <li>the {@link org.orekit.orbits.KeplerianOrbit Keplerian orbit parameters} (a, e, i, ω, Ω,
103  *   M or E or v) in meters and radians,</li>
104  *   <li>the {@link org.orekit.orbits.CircularOrbit circular orbit parameters} (a, e<sub>x</sub>, e<sub>y</sub>, i,
105  *   Ω, α<sub>M</sub> or α<sub>E</sub> or α<sub>v</sub>) in meters
106  *   and radians,</li>
107  *   <li>the {@link org.orekit.orbits.CartesianOrbit Cartesian orbit parameters} (x, y, z, v<sub>x</sub>,
108  *   v<sub>y</sub>, v<sub>z</sub>) in meters and meters per seconds.
109  * </ul>
110  * <p> The last element is the mass in kilograms.
111  *
112  * <p>The following code snippet shows a typical setting for Low Earth Orbit propagation in
113  * equinoctial parameters and true longitude argument:</p>
114  * <pre>
115  * final double dP       = 0.001;
116  * final double minStep  = 0.001;
117  * final double maxStep  = 500;
118  * final double initStep = 60;
119  * final double[][] tolerance = NumericalPropagator.tolerances(dP, orbit, OrbitType.EQUINOCTIAL);
120  * AdaptiveStepsizeIntegrator integrator = new DormandPrince853Integrator(minStep, maxStep, tolerance[0], tolerance[1]);
121  * integrator.setInitialStepSize(initStep);
122  * propagator = new NumericalPropagator(integrator);
123  * </pre>
124  * <p>By default, at the end of the propagation, the propagator resets the initial state to the final state,
125  * thus allowing a new propagation to be started from there without recomputing the part already performed.
126  * This behaviour can be chenged by calling {@link #setResetAtEnd(boolean)}.
127  * </p>
128  * <p>Beware the same instance cannot be used simultaneously by different threads, the class is <em>not</em>
129  * thread-safe.</p>
130  *
131  * @see SpacecraftState
132  * @see ForceModel
133  * @see org.orekit.propagation.sampling.OrekitStepHandler
134  * @see org.orekit.propagation.sampling.OrekitFixedStepHandler
135  * @see org.orekit.propagation.integration.IntegratedEphemeris
136  * @see TimeDerivativesEquations
137  *
138  * @author Mathieu Rom&eacute;ro
139  * @author Luc Maisonobe
140  * @author Guylaine Prat
141  * @author Fabien Maussion
142  * @author V&eacute;ronique Pommier-Maurussane
143  */
144 public class NumericalPropagator extends AbstractIntegratedPropagator {
145 
146     /** Force models used during the extrapolation of the orbit. */
147     private final List<ForceModel> forceModels;
148 
149     /** boolean to ignore or not the creation of a NewtonianAttraction. */
150     private boolean ignoreCentralAttraction = false;
151 
152     /** Create a new instance of NumericalPropagator, based on orbit definition mu.
153      * After creation, the instance is empty, i.e. the attitude provider is set to an
154      * unspecified default law and there are no perturbing forces at all.
155      * This means that if {@link #addForceModel addForceModel} is not
156      * called after creation, the integrated orbit will follow a Keplerian
157      * evolution only. The defaults are {@link OrbitType#EQUINOCTIAL}
158      * for {@link #setOrbitType(OrbitType) propagation
159      * orbit type} and {@link PositionAngle#TRUE} for {@link
160      * #setPositionAngleType(PositionAngle) position angle type}.
161      *
162      * <p>This constructor uses the {@link DataContext#getDefault() default data context}.
163      *
164      * @param integrator numerical integrator to use for propagation.
165      * @see #NumericalPropagator(ODEIntegrator, AttitudeProvider)
166      */
167     @DefaultDataContext
168     public NumericalPropagator(final ODEIntegrator integrator) {
169         this(integrator,
170                 Propagator.getDefaultLaw(DataContext.getDefault().getFrames()));
171     }
172 
173     /** Create a new instance of NumericalPropagator, based on orbit definition mu.
174      * After creation, the instance is empty, i.e. the attitude provider is set to an
175      * unspecified default law and there are no perturbing forces at all.
176      * This means that if {@link #addForceModel addForceModel} is not
177      * called after creation, the integrated orbit will follow a Keplerian
178      * evolution only. The defaults are {@link OrbitType#EQUINOCTIAL}
179      * for {@link #setOrbitType(OrbitType) propagation
180      * orbit type} and {@link PositionAngle#TRUE} for {@link
181      * #setPositionAngleType(PositionAngle) position angle type}.
182      * @param integrator numerical integrator to use for propagation.
183      * @param attitudeProvider the attitude law.
184      * @since 10.1
185      */
186     public NumericalPropagator(final ODEIntegrator integrator,
187                                final AttitudeProvider attitudeProvider) {
188         super(integrator, PropagationType.MEAN);
189         forceModels = new ArrayList<ForceModel>();
190         initMapper();
191         setAttitudeProvider(attitudeProvider);
192         clearStepHandlers();
193         setOrbitType(OrbitType.EQUINOCTIAL);
194         setPositionAngleType(PositionAngle.TRUE);
195     }
196 
197     /** Set the flag to ignore or not the creation of a {@link NewtonianAttraction}.
198      * @param ignoreCentralAttraction if true, {@link NewtonianAttraction} is <em>not</em>
199      * added automatically if missing
200      */
201     public void setIgnoreCentralAttraction(final boolean ignoreCentralAttraction) {
202         this.ignoreCentralAttraction = ignoreCentralAttraction;
203     }
204 
205      /** Set the central attraction coefficient μ.
206       * <p>
207       * Setting the central attraction coefficient is
208       * equivalent to {@link #addForceModel(ForceModel) add}
209       * a {@link NewtonianAttraction} force model.
210       * </p>
211      * @param mu central attraction coefficient (m³/s²)
212      * @see #addForceModel(ForceModel)
213      * @see #getAllForceModels()
214      */
215     public void setMu(final double mu) {
216         if (ignoreCentralAttraction) {
217             superSetMu(mu);
218         } else {
219             addForceModel(new NewtonianAttraction(mu));
220         }
221     }
222 
223     /** Set the central attraction coefficient μ only in upper class.
224      * @param mu central attraction coefficient (m³/s²)
225      */
226     private void superSetMu(final double mu) {
227         super.setMu(mu);
228     }
229 
230     /** Check if Newtonian attraction force model is available.
231      * <p>
232      * Newtonian attraction is always the last force model in the list.
233      * </p>
234      * @return true if Newtonian attraction force model is available
235      */
236     private boolean hasNewtonianAttraction() {
237         final int last = forceModels.size() - 1;
238         return last >= 0 && forceModels.get(last) instanceof NewtonianAttraction;
239     }
240 
241     /** Add a force model.
242      * <p>If this method is not called at all, the integrated orbit will follow
243      * a Keplerian evolution only.</p>
244      * @param model {@link ForceModel} to add (it can be either a perturbing force
245      * model or an instance of {@link NewtonianAttraction})
246      * @see #removeForceModels()
247      * @see #setMu(double)
248      */
249     public void addForceModel(final ForceModel model) {
250 
251         if (model instanceof NewtonianAttraction) {
252             // we want to add the central attraction force model
253 
254             try {
255                 // ensure we are notified of any mu change
256                 model.getParametersDrivers().get(0).addObserver(new ParameterObserver() {
257                     /** {@inheritDoc} */
258                     @Override
259                     public void valueChanged(final double previousValue, final ParameterDriver driver) {
260                         superSetMu(driver.getValue());
261                     }
262                 });
263             } catch (OrekitException oe) {
264                 // this should never happen
265                 throw new OrekitInternalError(oe);
266             }
267 
268             if (hasNewtonianAttraction()) {
269                 // there is already a central attraction model, replace it
270                 forceModels.set(forceModels.size() - 1, model);
271             } else {
272                 // there are no central attraction model yet, add it at the end of the list
273                 forceModels.add(model);
274             }
275         } else {
276             // we want to add a perturbing force model
277             if (hasNewtonianAttraction()) {
278                 // insert the new force model before Newtonian attraction,
279                 // which should always be the last one in the list
280                 forceModels.add(forceModels.size() - 1, model);
281             } else {
282                 // we only have perturbing force models up to now, just append at the end of the list
283                 forceModels.add(model);
284             }
285         }
286 
287     }
288 
289     /** Remove all force models (except central attraction).
290      * <p>Once all perturbing forces have been removed (and as long as no new force
291      * model is added), the integrated orbit will follow a Keplerian evolution
292      * only.</p>
293      * @see #addForceModel(ForceModel)
294      */
295     public void removeForceModels() {
296         final int last = forceModels.size() - 1;
297         if (hasNewtonianAttraction()) {
298             // preserve the Newtonian attraction model at the end
299             final ForceModel newton = forceModels.get(last);
300             forceModels.clear();
301             forceModels.add(newton);
302         } else {
303             forceModels.clear();
304         }
305     }
306 
307     /** Get all the force models, perturbing forces and Newtonian attraction included.
308      * @return list of perturbing force models, with Newtonian attraction being the
309      * last one
310      * @see #addForceModel(ForceModel)
311      * @see #setMu(double)
312      */
313     public List<ForceModel> getAllForceModels() {
314         return Collections.unmodifiableList(forceModels);
315     }
316 
317     /** Set propagation orbit type.
318      * @param orbitType orbit type to use for propagation, null for
319      * propagating using {@link org.orekit.utils.AbsolutePVCoordinates} rather than {@link Orbit}
320      */
321     public void setOrbitType(final OrbitType orbitType) {
322         super.setOrbitType(orbitType);
323     }
324 
325     /** Get propagation parameter type.
326      * @return orbit type used for propagation, null for
327      * propagating using {@link org.orekit.utils.AbsolutePVCoordinates} rather than {@link Orbit}
328      */
329     public OrbitType getOrbitType() {
330         return super.getOrbitType();
331     }
332 
333     /** Set position angle type.
334      * <p>
335      * The position parameter type is meaningful only if {@link
336      * #getOrbitType() propagation orbit type}
337      * support it. As an example, it is not meaningful for propagation
338      * in {@link OrbitType#CARTESIAN Cartesian} parameters.
339      * </p>
340      * @param positionAngleType angle type to use for propagation
341      */
342     public void setPositionAngleType(final PositionAngle positionAngleType) {
343         super.setPositionAngleType(positionAngleType);
344     }
345 
346     /** Get propagation parameter type.
347      * @return angle type to use for propagation
348      */
349     public PositionAngle getPositionAngleType() {
350         return super.getPositionAngleType();
351     }
352 
353     /** Set the initial state.
354      * @param initialState initial state
355      */
356     public void setInitialState(final SpacecraftState initialState) {
357         resetInitialState(initialState);
358     }
359 
360     /** {@inheritDoc} */
361     public void resetInitialState(final SpacecraftState state) {
362         super.resetInitialState(state);
363         if (!hasNewtonianAttraction()) {
364             // use the state to define central attraction
365             setMu(state.getMu());
366         }
367         setStartDate(state.getDate());
368     }
369 
370     /** {@inheritDoc} */
371     public TimeStampedPVCoordinates getPVCoordinates(final AbsoluteDate date, final Frame frame) {
372         return propagate(date).getPVCoordinates(frame);
373     }
374 
375     /** {@inheritDoc} */
376     protected StateMapper createMapper(final AbsoluteDate referenceDate, final double mu,
377                                        final OrbitType orbitType, final PositionAngle positionAngleType,
378                                        final AttitudeProvider attitudeProvider, final Frame frame) {
379         return new OsculatingMapper(referenceDate, mu, orbitType, positionAngleType, attitudeProvider, frame);
380     }
381 
382     /** Internal mapper using directly osculating parameters. */
383     private static class OsculatingMapper extends StateMapper {
384 
385         /** Simple constructor.
386          * <p>
387          * The position parameter type is meaningful only if {@link
388          * #getOrbitType() propagation orbit type}
389          * support it. As an example, it is not meaningful for propagation
390          * in {@link OrbitType#CARTESIAN Cartesian} parameters.
391          * </p>
392          * @param referenceDate reference date
393          * @param mu central attraction coefficient (m³/s²)
394          * @param orbitType orbit type to use for mapping (can be null for {@link AbsolutePVCoordinates})
395          * @param positionAngleType angle type to use for propagation
396          * @param attitudeProvider attitude provider
397          * @param frame inertial frame
398          */
399         OsculatingMapper(final AbsoluteDate referenceDate, final double mu,
400                          final OrbitType orbitType, final PositionAngle positionAngleType,
401                          final AttitudeProvider attitudeProvider, final Frame frame) {
402             super(referenceDate, mu, orbitType, positionAngleType, attitudeProvider, frame);
403         }
404 
405         /** {@inheritDoc} */
406         public SpacecraftState mapArrayToState(final AbsoluteDate date, final double[] y, final double[] yDot,
407                                                final PropagationType type) {
408             // the parameter type is ignored for the Numerical Propagator
409 
410             final double mass = y[6];
411             if (mass <= 0.0) {
412                 throw new OrekitException(OrekitMessages.SPACECRAFT_MASS_BECOMES_NEGATIVE, mass);
413             }
414 
415             if (getOrbitType() == null) {
416                 // propagation uses absolute position-velocity-acceleration
417                 final Vector3D p = new Vector3D(y[0],    y[1],    y[2]);
418                 final Vector3D v = new Vector3D(y[3],    y[4],    y[5]);
419                 final Vector3D a;
420                 final AbsolutePVCoordinates absPva;
421                 if (yDot == null) {
422                     absPva = new AbsolutePVCoordinates(getFrame(), new TimeStampedPVCoordinates(date, p, v));
423                 } else {
424                     a = new Vector3D(yDot[3], yDot[4], yDot[5]);
425                     absPva = new AbsolutePVCoordinates(getFrame(), new TimeStampedPVCoordinates(date, p, v, a));
426                 }
427 
428                 final Attitude attitude = getAttitudeProvider().getAttitude(absPva, date, getFrame());
429                 return new SpacecraftState(absPva, attitude, mass);
430             } else {
431                 // propagation uses regular orbits
432                 final Orbit orbit       = getOrbitType().mapArrayToOrbit(y, yDot, getPositionAngleType(), date, getMu(), getFrame());
433                 final Attitude attitude = getAttitudeProvider().getAttitude(orbit, date, getFrame());
434 
435                 return new SpacecraftState(orbit, attitude, mass);
436             }
437 
438         }
439 
440         /** {@inheritDoc} */
441         public void mapStateToArray(final SpacecraftState state, final double[] y, final double[] yDot) {
442             if (getOrbitType() == null) {
443                 // propagation uses absolute position-velocity-acceleration
444                 final Vector3D p = state.getAbsPVA().getPosition();
445                 final Vector3D v = state.getAbsPVA().getVelocity();
446                 y[0] = p.getX();
447                 y[1] = p.getY();
448                 y[2] = p.getZ();
449                 y[3] = v.getX();
450                 y[4] = v.getY();
451                 y[5] = v.getZ();
452                 y[6] = state.getMass();
453             }
454             else {
455                 getOrbitType().mapOrbitToArray(state.getOrbit(), getPositionAngleType(), y, yDot);
456                 y[6] = state.getMass();
457             }
458         }
459 
460     }
461 
462     /** {@inheritDoc} */
463     protected MainStateEquations getMainStateEquations(final ODEIntegrator integrator) {
464         return new Main(integrator);
465     }
466 
467     /** Internal class for osculating parameters integration. */
468     private class Main implements MainStateEquations, TimeDerivativesEquations {
469 
470         /** Derivatives array. */
471         private final double[] yDot;
472 
473         /** Current state. */
474         private SpacecraftState currentState;
475 
476         /** Jacobian of the orbital parameters with respect to the Cartesian parameters. */
477         private double[][] jacobian;
478 
479         /** Simple constructor.
480          * @param integrator numerical integrator to use for propagation.
481          */
482         Main(final ODEIntegrator integrator) {
483 
484             this.yDot     = new double[7];
485             this.jacobian = new double[6][6];
486 
487             for (final ForceModel forceModel : forceModels) {
488                 forceModel.getEventsDetectors().forEach(detector -> setUpEventDetector(integrator, detector));
489             }
490 
491             if (getOrbitType() == null) {
492                 // propagation uses absolute position-velocity-acceleration
493                 // we can set Jacobian once and for all
494                 for (int i = 0; i < jacobian.length; ++i) {
495                     Arrays.fill(jacobian[i], 0.0);
496                     jacobian[i][i] = 1.0;
497                 }
498             }
499 
500         }
501 
502         /** {@inheritDoc} */
503         @Override
504         public void init(final SpacecraftState initialState, final AbsoluteDate target) {
505             for (final ForceModel forceModel : forceModels) {
506                 forceModel.init(initialState, target);
507             }
508         }
509 
510         /** {@inheritDoc} */
511         @Override
512         public double[] computeDerivatives(final SpacecraftState state) {
513 
514             currentState = state;
515             Arrays.fill(yDot, 0.0);
516             if (getOrbitType() != null) {
517                 // propagation uses regular orbits
518                 currentState.getOrbit().getJacobianWrtCartesian(getPositionAngleType(), jacobian);
519             }
520 
521             // compute the contributions of all perturbing forces,
522             // using the Kepler contribution at the end since
523             // NewtonianAttraction is always the last instance in the list
524             for (final ForceModel forceModel : forceModels) {
525                 forceModel.addContribution(state, this);
526             }
527 
528             if (getOrbitType() == null) {
529                 // position derivative is velocity, and was not added above in the force models
530                 // (it is added when orbit type is non-null because NewtonianAttraction considers it)
531                 final Vector3D velocity = currentState.getPVCoordinates().getVelocity();
532                 yDot[0] += velocity.getX();
533                 yDot[1] += velocity.getY();
534                 yDot[2] += velocity.getZ();
535             }
536 
537 
538             return yDot.clone();
539 
540         }
541 
542         /** {@inheritDoc} */
543         @Override
544         public void addKeplerContribution(final double mu) {
545             if (getOrbitType() == null) {
546 
547                 // if mu is neither 0 nor NaN, we want to include Newtonian acceleration
548                 if (mu > 0) {
549                     // velocity derivative is Newtonian acceleration
550                     final Vector3D position = currentState.getPVCoordinates().getPosition();
551                     final double r2         = position.getNormSq();
552                     final double coeff      = -mu / (r2 * FastMath.sqrt(r2));
553                     yDot[3] += coeff * position.getX();
554                     yDot[4] += coeff * position.getY();
555                     yDot[5] += coeff * position.getZ();
556                 }
557 
558             } else {
559                 // propagation uses regular orbits
560                 currentState.getOrbit().addKeplerContribution(getPositionAngleType(), mu, yDot);
561             }
562         }
563 
564         /** {@inheritDoc} */
565         public void addNonKeplerianAcceleration(final Vector3D gamma) {
566             for (int i = 0; i < 6; ++i) {
567                 final double[] jRow = jacobian[i];
568                 yDot[i] += jRow[3] * gamma.getX() + jRow[4] * gamma.getY() + jRow[5] * gamma.getZ();
569             }
570         }
571 
572         /** {@inheritDoc} */
573         @Override
574         public void addMassDerivative(final double q) {
575             if (q > 0) {
576                 throw new OrekitIllegalArgumentException(OrekitMessages.POSITIVE_FLOW_RATE, q);
577             }
578             yDot[6] += q;
579         }
580 
581     }
582 
583     /** Estimate tolerance vectors for integrators when propagating in absolute position-velocity-acceleration.
584      * @param dP user specified position error
585      * @param absPva reference absolute position-velocity-acceleration
586      * @return a two rows array, row 0 being the absolute tolerance error and row 1
587      * being the relative tolerance error
588      * @see NumericalPropagator#tolerances(double, Orbit, OrbitType)
589      */
590     public static double[][] tolerances(final double dP, final AbsolutePVCoordinates absPva) {
591 
592         final double relative = dP / absPva.getPosition().getNorm();
593         final double dV = relative * absPva.getVelocity().getNorm();
594 
595         final double[] absTol = new double[7];
596         final double[] relTol = new double[7];
597 
598         absTol[0] = dP;
599         absTol[1] = dP;
600         absTol[2] = dP;
601         absTol[3] = dV;
602         absTol[4] = dV;
603         absTol[5] = dV;
604 
605         // we set the mass tolerance arbitrarily to 1.0e-6 kg, as mass evolves linearly
606         // with trust, this often has no influence at all on propagation
607         absTol[6] = 1.0e-6;
608 
609         Arrays.fill(relTol, relative);
610 
611         return new double[][] {
612             absTol, relTol
613         };
614 
615     }
616 
617     /** Estimate tolerance vectors for integrators when propagating in orbits.
618      * <p>
619      * The errors are estimated from partial derivatives properties of orbits,
620      * starting from a scalar position error specified by the user.
621      * Considering the energy conservation equation V = sqrt(mu (2/r - 1/a)),
622      * we get at constant energy (i.e. on a Keplerian trajectory):
623      * <pre>
624      * V² r |dV| = mu |dr|
625      * </pre>
626      * <p> So we deduce a scalar velocity error consistent with the position error.
627      * From here, we apply orbits Jacobians matrices to get consistent errors
628      * on orbital parameters.
629      *
630      * <p>
631      * The tolerances are only <em>orders of magnitude</em>, and integrator tolerances
632      * are only local estimates, not global ones. So some care must be taken when using
633      * these tolerances. Setting 1mm as a position error does NOT mean the tolerances
634      * will guarantee a 1mm error position after several orbits integration.
635      * </p>
636      * @param dP user specified position error
637      * @param orbit reference orbit
638      * @param type propagation type for the meaning of the tolerance vectors elements
639      * (it may be different from {@code orbit.getType()})
640      * @return a two rows array, row 0 being the absolute tolerance error and row 1
641      * being the relative tolerance error
642      */
643     public static double[][] tolerances(final double dP, final Orbit orbit, final OrbitType type) {
644 
645         // estimate the scalar velocity error
646         final PVCoordinates pv = orbit.getPVCoordinates();
647         final double r2 = pv.getPosition().getNormSq();
648         final double v  = pv.getVelocity().getNorm();
649         final double dV = orbit.getMu() * dP / (v * r2);
650 
651         return tolerances(dP, dV, orbit, type);
652 
653     }
654 
655     /** Estimate tolerance vectors for integrators when propagating in orbits.
656      * <p>
657      * The errors are estimated from partial derivatives properties of orbits,
658      * starting from scalar position and velocity errors specified by the user.
659      * <p>
660      * The tolerances are only <em>orders of magnitude</em>, and integrator tolerances
661      * are only local estimates, not global ones. So some care must be taken when using
662      * these tolerances. Setting 1mm as a position error does NOT mean the tolerances
663      * will guarantee a 1mm error position after several orbits integration.
664      * </p>
665      * @param dP user specified position error
666      * @param dV user specified velocity error
667      * @param orbit reference orbit
668      * @param type propagation type for the meaning of the tolerance vectors elements
669      * (it may be different from {@code orbit.getType()})
670      * @return a two rows array, row 0 being the absolute tolerance error and row 1
671      * being the relative tolerance error
672      * @since 10.3
673      */
674     public static double[][] tolerances(final double dP, final double dV,
675                                         final Orbit orbit, final OrbitType type) {
676 
677         final double[] absTol = new double[7];
678         final double[] relTol = new double[7];
679 
680         // we set the mass tolerance arbitrarily to 1.0e-6 kg, as mass evolves linearly
681         // with trust, this often has no influence at all on propagation
682         absTol[6] = 1.0e-6;
683 
684         if (type == OrbitType.CARTESIAN) {
685             absTol[0] = dP;
686             absTol[1] = dP;
687             absTol[2] = dP;
688             absTol[3] = dV;
689             absTol[4] = dV;
690             absTol[5] = dV;
691         } else {
692 
693             // convert the orbit to the desired type
694             final double[][] jacobian = new double[6][6];
695             final Orbit converted = type.convertType(orbit);
696             converted.getJacobianWrtCartesian(PositionAngle.TRUE, jacobian);
697 
698             for (int i = 0; i < 6; ++i) {
699                 final double[] row = jacobian[i];
700                 absTol[i] = FastMath.abs(row[0]) * dP +
701                             FastMath.abs(row[1]) * dP +
702                             FastMath.abs(row[2]) * dP +
703                             FastMath.abs(row[3]) * dV +
704                             FastMath.abs(row[4]) * dV +
705                             FastMath.abs(row[5]) * dV;
706                 if (Double.isNaN(absTol[i])) {
707                     throw new OrekitException(OrekitMessages.SINGULAR_JACOBIAN_FOR_ORBIT_TYPE, type);
708                 }
709             }
710 
711         }
712 
713         Arrays.fill(relTol, dP / FastMath.sqrt(orbit.getPVCoordinates().getPosition().getNormSq()));
714 
715         return new double[][] {
716             absTol, relTol
717         };
718 
719     }
720 
721     /** {@inheritDoc} */
722     @Override
723     protected void beforeIntegration(final SpacecraftState initialState, final AbsoluteDate tEnd) {
724 
725         if (!getFrame().isPseudoInertial()) {
726 
727             // inspect all force models to find InertialForces
728             for (ForceModel force : forceModels) {
729                 if (force instanceof InertialForces) {
730                     return;
731                 }
732             }
733 
734             // throw exception if no inertial forces found
735             throw new OrekitException(OrekitMessages.INERTIAL_FORCE_MODEL_MISSING, getFrame().getName());
736 
737         }
738 
739     }
740 
741 }