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