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