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