NumericalPropagator.java

  1. /* Copyright 2002-2020 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. import java.util.ArrayList;
  19. import java.util.Arrays;
  20. import java.util.Collections;
  21. import java.util.List;

  22. import org.hipparchus.geometry.euclidean.threed.Vector3D;
  23. import org.hipparchus.ode.ODEIntegrator;
  24. import org.hipparchus.util.FastMath;
  25. import org.orekit.annotation.DefaultDataContext;
  26. import org.orekit.attitudes.Attitude;
  27. import org.orekit.attitudes.AttitudeProvider;
  28. import org.orekit.data.DataContext;
  29. import org.orekit.errors.OrekitException;
  30. import org.orekit.errors.OrekitIllegalArgumentException;
  31. import org.orekit.errors.OrekitInternalError;
  32. import org.orekit.errors.OrekitMessages;
  33. import org.orekit.forces.ForceModel;
  34. import org.orekit.forces.gravity.NewtonianAttraction;
  35. import org.orekit.forces.inertia.InertialForces;
  36. import org.orekit.frames.Frame;
  37. import org.orekit.orbits.Orbit;
  38. import org.orekit.orbits.OrbitType;
  39. import org.orekit.orbits.PositionAngle;
  40. import org.orekit.propagation.PropagationType;
  41. import org.orekit.propagation.Propagator;
  42. import org.orekit.propagation.SpacecraftState;
  43. import org.orekit.propagation.events.EventDetector;
  44. import org.orekit.propagation.integration.AbstractIntegratedPropagator;
  45. import org.orekit.propagation.integration.StateMapper;
  46. import org.orekit.time.AbsoluteDate;
  47. import org.orekit.utils.AbsolutePVCoordinates;
  48. import org.orekit.utils.PVCoordinates;
  49. import org.orekit.utils.ParameterDriver;
  50. import org.orekit.utils.ParameterObserver;
  51. import org.orekit.utils.TimeStampedPVCoordinates;

  52. /** This class propagates {@link org.orekit.orbits.Orbit orbits} using
  53.  * numerical integration.
  54.  * <p>Numerical propagation is much more accurate than analytical propagation
  55.  * like for example {@link org.orekit.propagation.analytical.KeplerianPropagator
  56.  * Keplerian} or {@link org.orekit.propagation.analytical.EcksteinHechlerPropagator
  57.  * Eckstein-Hechler}, but requires a few more steps to set up to be used properly.
  58.  * Whereas analytical propagators are configured only thanks to their various
  59.  * constructors and can be used immediately after construction, numerical propagators
  60.  * configuration involve setting several parameters between construction time
  61.  * and propagation time.</p>
  62.  * <p>The configuration parameters that can be set are:</p>
  63.  * <ul>
  64.  *   <li>the initial spacecraft state ({@link #setInitialState(SpacecraftState)})</li>
  65.  *   <li>the central attraction coefficient ({@link #setMu(double)})</li>
  66.  *   <li>the various force models ({@link #addForceModel(ForceModel)},
  67.  *   {@link #removeForceModels()})</li>
  68.  *   <li>the {@link OrbitType type} of orbital parameters to be used for propagation
  69.  *   ({@link #setOrbitType(OrbitType)}),
  70.  *   <li>the {@link PositionAngle type} of position angle to be used in orbital parameters
  71.  *   to be used for propagation where it is relevant ({@link
  72.  *   #setPositionAngleType(PositionAngle)}),
  73.  *   <li>whether {@link org.orekit.propagation.integration.AdditionalEquations additional equations}
  74.  *   (for example {@link PartialDerivativesEquations Jacobians}) should be propagated along with orbital state
  75.  *   ({@link #addAdditionalEquations(org.orekit.propagation.integration.AdditionalEquations)}),
  76.  *   <li>the discrete events that should be triggered during propagation
  77.  *   ({@link #addEventDetector(EventDetector)},
  78.  *   {@link #clearEventsDetectors()})</li>
  79.  *   <li>the binding logic with the rest of the application ({@link #setSlaveMode()},
  80.  *   {@link #setMasterMode(double, org.orekit.propagation.sampling.OrekitFixedStepHandler)},
  81.  *   {@link #setMasterMode(org.orekit.propagation.sampling.OrekitStepHandler)},
  82.  *   {@link #setEphemerisMode()}, {@link #getGeneratedEphemeris()})</li>
  83.  * </ul>
  84.  * <p>From these configuration parameters, only the initial state is mandatory. The default
  85.  * propagation settings are in {@link OrbitType#EQUINOCTIAL equinoctial} parameters with
  86.  * {@link PositionAngle#TRUE true} longitude argument. If the central attraction coefficient
  87.  * is not explicitly specified, the one used to define the initial orbit will be used.
  88.  * However, specifying only the initial state and perhaps the central attraction coefficient
  89.  * would mean the propagator would use only Keplerian forces. In this case, the simpler {@link
  90.  * org.orekit.propagation.analytical.KeplerianPropagator KeplerianPropagator} class would
  91.  * perhaps be more effective.</p>
  92.  * <p>The underlying numerical integrator set up in the constructor may also have its own
  93.  * configuration parameters. Typical configuration parameters for adaptive stepsize integrators
  94.  * are the min, max and perhaps start step size as well as the absolute and/or relative errors
  95.  * thresholds.</p>
  96.  * <p>The state that is seen by the integrator is a simple seven elements double array.
  97.  * The six first elements are either:
  98.  * <ul>
  99.  *   <li>the {@link org.orekit.orbits.EquinoctialOrbit equinoctial orbit parameters} (a, e<sub>x</sub>,
  100.  *   e<sub>y</sub>, h<sub>x</sub>, h<sub>y</sub>, λ<sub>M</sub> or λ<sub>E</sub>
  101.  *   or λ<sub>v</sub>) in meters and radians,</li>
  102.  *   <li>the {@link org.orekit.orbits.KeplerianOrbit Keplerian orbit parameters} (a, e, i, ω, Ω,
  103.  *   M or E or v) in meters and radians,</li>
  104.  *   <li>the {@link org.orekit.orbits.CircularOrbit circular orbit parameters} (a, e<sub>x</sub>, e<sub>y</sub>, i,
  105.  *   Ω, α<sub>M</sub> or α<sub>E</sub> or α<sub>v</sub>) in meters
  106.  *   and radians,</li>
  107.  *   <li>the {@link org.orekit.orbits.CartesianOrbit Cartesian orbit parameters} (x, y, z, v<sub>x</sub>,
  108.  *   v<sub>y</sub>, v<sub>z</sub>) in meters and meters per seconds.
  109.  * </ul>
  110.  * <p> The last element is the mass in kilograms.
  111.  *
  112.  * <p>The following code snippet shows a typical setting for Low Earth Orbit propagation in
  113.  * equinoctial parameters and true longitude argument:</p>
  114.  * <pre>
  115.  * final double dP       = 0.001;
  116.  * final double minStep  = 0.001;
  117.  * final double maxStep  = 500;
  118.  * final double initStep = 60;
  119.  * final double[][] tolerance = NumericalPropagator.tolerances(dP, orbit, OrbitType.EQUINOCTIAL);
  120.  * AdaptiveStepsizeIntegrator integrator = new DormandPrince853Integrator(minStep, maxStep, tolerance[0], tolerance[1]);
  121.  * integrator.setInitialStepSize(initStep);
  122.  * propagator = new NumericalPropagator(integrator);
  123.  * </pre>
  124.  * <p>By default, at the end of the propagation, the propagator resets the initial state to the final state,
  125.  * thus allowing a new propagation to be started from there without recomputing the part already performed.
  126.  * This behaviour can be chenged by calling {@link #setResetAtEnd(boolean)}.
  127.  * </p>
  128.  * <p>Beware the same instance cannot be used simultaneously by different threads, the class is <em>not</em>
  129.  * thread-safe.</p>
  130.  *
  131.  * @see SpacecraftState
  132.  * @see ForceModel
  133.  * @see org.orekit.propagation.sampling.OrekitStepHandler
  134.  * @see org.orekit.propagation.sampling.OrekitFixedStepHandler
  135.  * @see org.orekit.propagation.integration.IntegratedEphemeris
  136.  * @see TimeDerivativesEquations
  137.  *
  138.  * @author Mathieu Rom&eacute;ro
  139.  * @author Luc Maisonobe
  140.  * @author Guylaine Prat
  141.  * @author Fabien Maussion
  142.  * @author V&eacute;ronique Pommier-Maurussane
  143.  */
  144. public class NumericalPropagator extends AbstractIntegratedPropagator {

  145.     /** Force models used during the extrapolation of the orbit. */
  146.     private final List<ForceModel> forceModels;

  147.     /** boolean to ignore or not the creation of a NewtonianAttraction. */
  148.     private boolean ignoreCentralAttraction = false;

  149.     /** Create a new instance of NumericalPropagator, based on orbit definition mu.
  150.      * After creation, the instance is empty, i.e. the attitude provider is set to an
  151.      * unspecified default law and there are no perturbing forces at all.
  152.      * This means that if {@link #addForceModel addForceModel} is not
  153.      * called after creation, the integrated orbit will follow a Keplerian
  154.      * evolution only. The defaults are {@link OrbitType#EQUINOCTIAL}
  155.      * for {@link #setOrbitType(OrbitType) propagation
  156.      * orbit type} and {@link PositionAngle#TRUE} for {@link
  157.      * #setPositionAngleType(PositionAngle) position angle type}.
  158.      *
  159.      * <p>This constructor uses the {@link DataContext#getDefault() default data context}.
  160.      *
  161.      * @param integrator numerical integrator to use for propagation.
  162.      * @see #NumericalPropagator(ODEIntegrator, AttitudeProvider)
  163.      */
  164.     @DefaultDataContext
  165.     public NumericalPropagator(final ODEIntegrator integrator) {
  166.         this(integrator,
  167.                 Propagator.getDefaultLaw(DataContext.getDefault().getFrames()));
  168.     }

  169.     /** Create a new instance of NumericalPropagator, based on orbit definition mu.
  170.      * After creation, the instance is empty, i.e. the attitude provider is set to an
  171.      * unspecified default law and there are no perturbing forces at all.
  172.      * This means that if {@link #addForceModel addForceModel} is not
  173.      * called after creation, the integrated orbit will follow a Keplerian
  174.      * evolution only. The defaults are {@link OrbitType#EQUINOCTIAL}
  175.      * for {@link #setOrbitType(OrbitType) propagation
  176.      * orbit type} and {@link PositionAngle#TRUE} for {@link
  177.      * #setPositionAngleType(PositionAngle) position angle type}.
  178.      * @param integrator numerical integrator to use for propagation.
  179.      * @param attitudeProvider the attitude law.
  180.      * @since 10.1
  181.      */
  182.     public NumericalPropagator(final ODEIntegrator integrator,
  183.                                final AttitudeProvider attitudeProvider) {
  184.         super(integrator, PropagationType.MEAN);
  185.         forceModels = new ArrayList<ForceModel>();
  186.         initMapper();
  187.         setAttitudeProvider(attitudeProvider);
  188.         setSlaveMode();
  189.         setOrbitType(OrbitType.EQUINOCTIAL);
  190.         setPositionAngleType(PositionAngle.TRUE);
  191.     }

  192.     /** Set the flag to ignore or not the creation of a {@link NewtonianAttraction}.
  193.      * @param ignoreCentralAttraction if true, {@link NewtonianAttraction} is <em>not</em>
  194.      * added automatically if missing
  195.      */
  196.     public void setIgnoreCentralAttraction(final boolean ignoreCentralAttraction) {
  197.         this.ignoreCentralAttraction = ignoreCentralAttraction;
  198.     }

  199.      /** Set the central attraction coefficient μ.
  200.       * <p>
  201.       * Setting the central attraction coefficient is
  202.       * equivalent to {@link #addForceModel(ForceModel) add}
  203.       * a {@link NewtonianAttraction} force model.
  204.       * </p>
  205.      * @param mu central attraction coefficient (m³/s²)
  206.      * @see #addForceModel(ForceModel)
  207.      * @see #getAllForceModels()
  208.      */
  209.     public void setMu(final double mu) {
  210.         if (ignoreCentralAttraction) {
  211.             superSetMu(mu);
  212.         } else {
  213.             addForceModel(new NewtonianAttraction(mu));
  214.         }
  215.     }

  216.     /** Set the central attraction coefficient μ only in upper class.
  217.      * @param mu central attraction coefficient (m³/s²)
  218.      */
  219.     private void superSetMu(final double mu) {
  220.         super.setMu(mu);
  221.     }

  222.     /** Check if Newtonian attraction force model is available.
  223.      * <p>
  224.      * Newtonian attraction is always the last force model in the list.
  225.      * </p>
  226.      * @return true if Newtonian attraction force model is available
  227.      */
  228.     private boolean hasNewtonianAttraction() {
  229.         final int last = forceModels.size() - 1;
  230.         return last >= 0 && forceModels.get(last) instanceof NewtonianAttraction;
  231.     }

  232.     /** Add a force model.
  233.      * <p>If this method is not called at all, the integrated orbit will follow
  234.      * a Keplerian evolution only.</p>
  235.      * @param model {@link ForceModel} to add (it can be either a perturbing force
  236.      * model or an instance of {@link NewtonianAttraction})
  237.      * @see #removeForceModels()
  238.      * @see #setMu(double)
  239.      */
  240.     public void addForceModel(final ForceModel model) {

  241.         if (model instanceof NewtonianAttraction) {
  242.             // we want to add the central attraction force model

  243.             try {
  244.                 // ensure we are notified of any mu change
  245.                 model.getParametersDrivers()[0].addObserver(new ParameterObserver() {
  246.                     /** {@inheritDoc} */
  247.                     @Override
  248.                     public void valueChanged(final double previousValue, final ParameterDriver driver) {
  249.                         superSetMu(driver.getValue());
  250.                     }
  251.                 });
  252.             } catch (OrekitException oe) {
  253.                 // this should never happen
  254.                 throw new OrekitInternalError(oe);
  255.             }

  256.             if (hasNewtonianAttraction()) {
  257.                 // there is already a central attraction model, replace it
  258.                 forceModels.set(forceModels.size() - 1, model);
  259.             } else {
  260.                 // there are no central attraction model yet, add it at the end of the list
  261.                 forceModels.add(model);
  262.             }
  263.         } else {
  264.             // we want to add a perturbing force model
  265.             if (hasNewtonianAttraction()) {
  266.                 // insert the new force model before Newtonian attraction,
  267.                 // which should always be the last one in the list
  268.                 forceModels.add(forceModels.size() - 1, model);
  269.             } else {
  270.                 // we only have perturbing force models up to now, just append at the end of the list
  271.                 forceModels.add(model);
  272.             }
  273.         }

  274.     }

  275.     /** Remove all force models (except central attraction).
  276.      * <p>Once all perturbing forces have been removed (and as long as no new force
  277.      * model is added), the integrated orbit will follow a Keplerian evolution
  278.      * only.</p>
  279.      * @see #addForceModel(ForceModel)
  280.      */
  281.     public void removeForceModels() {
  282.         final int last = forceModels.size() - 1;
  283.         if (hasNewtonianAttraction()) {
  284.             // preserve the Newtonian attraction model at the end
  285.             final ForceModel newton = forceModels.get(last);
  286.             forceModels.clear();
  287.             forceModels.add(newton);
  288.         } else {
  289.             forceModels.clear();
  290.         }
  291.     }

  292.     /** Get all the force models, perturbing forces and Newtonian attraction included.
  293.      * @return list of perturbing force models, with Newtonian attraction being the
  294.      * last one
  295.      * @see #addForceModel(ForceModel)
  296.      * @see #setMu(double)
  297.      */
  298.     public List<ForceModel> getAllForceModels() {
  299.         return Collections.unmodifiableList(forceModels);
  300.     }

  301.     /** Set propagation orbit type.
  302.      * @param orbitType orbit type to use for propagation, null for
  303.      * propagating using {@link org.orekit.utils.AbsolutePVCoordinates} rather than {@link Orbit}
  304.      */
  305.     public void setOrbitType(final OrbitType orbitType) {
  306.         super.setOrbitType(orbitType);
  307.     }

  308.     /** Get propagation parameter type.
  309.      * @return orbit type used for propagation, null for
  310.      * propagating using {@link org.orekit.utils.AbsolutePVCoordinates} rather than {@link Orbit}
  311.      */
  312.     public OrbitType getOrbitType() {
  313.         return super.getOrbitType();
  314.     }

  315.     /** Set position angle type.
  316.      * <p>
  317.      * The position parameter type is meaningful only if {@link
  318.      * #getOrbitType() propagation orbit type}
  319.      * support it. As an example, it is not meaningful for propagation
  320.      * in {@link OrbitType#CARTESIAN Cartesian} parameters.
  321.      * </p>
  322.      * @param positionAngleType angle type to use for propagation
  323.      */
  324.     public void setPositionAngleType(final PositionAngle positionAngleType) {
  325.         super.setPositionAngleType(positionAngleType);
  326.     }

  327.     /** Get propagation parameter type.
  328.      * @return angle type to use for propagation
  329.      */
  330.     public PositionAngle getPositionAngleType() {
  331.         return super.getPositionAngleType();
  332.     }

  333.     /** Set the initial state.
  334.      * @param initialState initial state
  335.      */
  336.     public void setInitialState(final SpacecraftState initialState) {
  337.         resetInitialState(initialState);
  338.     }

  339.     /** {@inheritDoc} */
  340.     public void resetInitialState(final SpacecraftState state) {
  341.         super.resetInitialState(state);
  342.         if (!hasNewtonianAttraction()) {
  343.             // use the state to define central attraction
  344.             setMu(state.getMu());
  345.         }
  346.         setStartDate(state.getDate());
  347.     }

  348.     /** {@inheritDoc} */
  349.     public TimeStampedPVCoordinates getPVCoordinates(final AbsoluteDate date, final Frame frame) {
  350.         return propagate(date).getPVCoordinates(frame);
  351.     }

  352.     /** {@inheritDoc} */
  353.     protected StateMapper createMapper(final AbsoluteDate referenceDate, final double mu,
  354.                                        final OrbitType orbitType, final PositionAngle positionAngleType,
  355.                                        final AttitudeProvider attitudeProvider, final Frame frame) {
  356.         return new OsculatingMapper(referenceDate, mu, orbitType, positionAngleType, attitudeProvider, frame);
  357.     }

  358.     /** Internal mapper using directly osculating parameters. */
  359.     private static class OsculatingMapper extends StateMapper {

  360.         /** Simple constructor.
  361.          * <p>
  362.          * The position parameter type is meaningful only if {@link
  363.          * #getOrbitType() propagation orbit type}
  364.          * support it. As an example, it is not meaningful for propagation
  365.          * in {@link OrbitType#CARTESIAN Cartesian} parameters.
  366.          * </p>
  367.          * @param referenceDate reference date
  368.          * @param mu central attraction coefficient (m³/s²)
  369.          * @param orbitType orbit type to use for mapping (can be null for {@link AbsolutePVCoordinates})
  370.          * @param positionAngleType angle type to use for propagation
  371.          * @param attitudeProvider attitude provider
  372.          * @param frame inertial frame
  373.          */
  374.         OsculatingMapper(final AbsoluteDate referenceDate, final double mu,
  375.                          final OrbitType orbitType, final PositionAngle positionAngleType,
  376.                          final AttitudeProvider attitudeProvider, final Frame frame) {
  377.             super(referenceDate, mu, orbitType, positionAngleType, attitudeProvider, frame);
  378.         }

  379.         /** {@inheritDoc} */
  380.         public SpacecraftState mapArrayToState(final AbsoluteDate date, final double[] y, final double[] yDot,
  381.                                                final PropagationType type) {
  382.             // the parameter type is ignored for the Numerical Propagator

  383.             final double mass = y[6];
  384.             if (mass <= 0.0) {
  385.                 throw new OrekitException(OrekitMessages.SPACECRAFT_MASS_BECOMES_NEGATIVE, mass);
  386.             }

  387.             if (getOrbitType() == null) {
  388.                 // propagation uses absolute position-velocity-acceleration
  389.                 final Vector3D p = new Vector3D(y[0],    y[1],    y[2]);
  390.                 final Vector3D v = new Vector3D(y[3],    y[4],    y[5]);
  391.                 final Vector3D a;
  392.                 final AbsolutePVCoordinates absPva;
  393.                 if (yDot == null) {
  394.                     absPva = new AbsolutePVCoordinates(getFrame(), new TimeStampedPVCoordinates(date, p, v));
  395.                 } else {
  396.                     a = new Vector3D(yDot[3], yDot[4], yDot[5]);
  397.                     absPva = new AbsolutePVCoordinates(getFrame(), new TimeStampedPVCoordinates(date, p, v, a));
  398.                 }

  399.                 final Attitude attitude = getAttitudeProvider().getAttitude(absPva, date, getFrame());
  400.                 return new SpacecraftState(absPva, attitude, mass);
  401.             } else {
  402.                 // propagation uses regular orbits
  403.                 final Orbit orbit       = getOrbitType().mapArrayToOrbit(y, yDot, getPositionAngleType(), date, getMu(), getFrame());
  404.                 final Attitude attitude = getAttitudeProvider().getAttitude(orbit, date, getFrame());

  405.                 return new SpacecraftState(orbit, attitude, mass);
  406.             }

  407.         }

  408.         /** {@inheritDoc} */
  409.         public void mapStateToArray(final SpacecraftState state, final double[] y, final double[] yDot) {
  410.             if (getOrbitType() == null) {
  411.                 // propagation uses absolute position-velocity-acceleration
  412.                 final Vector3D p = state.getAbsPVA().getPosition();
  413.                 final Vector3D v = state.getAbsPVA().getVelocity();
  414.                 y[0] = p.getX();
  415.                 y[1] = p.getY();
  416.                 y[2] = p.getZ();
  417.                 y[3] = v.getX();
  418.                 y[4] = v.getY();
  419.                 y[5] = v.getZ();
  420.                 y[6] = state.getMass();
  421.             }
  422.             else {
  423.                 getOrbitType().mapOrbitToArray(state.getOrbit(), getPositionAngleType(), y, yDot);
  424.                 y[6] = state.getMass();
  425.             }
  426.         }

  427.     }

  428.     /** {@inheritDoc} */
  429.     protected MainStateEquations getMainStateEquations(final ODEIntegrator integrator) {
  430.         return new Main(integrator);
  431.     }

  432.     /** Internal class for osculating parameters integration. */
  433.     private class Main implements MainStateEquations, TimeDerivativesEquations {

  434.         /** Derivatives array. */
  435.         private final double[] yDot;

  436.         /** Current state. */
  437.         private SpacecraftState currentState;

  438.         /** Jacobian of the orbital parameters with respect to the Cartesian parameters. */
  439.         private double[][] jacobian;

  440.         /** Simple constructor.
  441.          * @param integrator numerical integrator to use for propagation.
  442.          */
  443.         Main(final ODEIntegrator integrator) {

  444.             this.yDot     = new double[7];
  445.             this.jacobian = new double[6][6];

  446.             for (final ForceModel forceModel : forceModels) {
  447.                 forceModel.getEventsDetectors().forEach(detector -> setUpEventDetector(integrator, detector));
  448.             }

  449.             if (getOrbitType() == null) {
  450.                 // propagation uses absolute position-velocity-acceleration
  451.                 // we can set Jacobian once and for all
  452.                 for (int i = 0; i < jacobian.length; ++i) {
  453.                     Arrays.fill(jacobian[i], 0.0);
  454.                     jacobian[i][i] = 1.0;
  455.                 }
  456.             }

  457.         }

  458.         /** {@inheritDoc} */
  459.         @Override
  460.         public void init(final SpacecraftState initialState, final AbsoluteDate target) {
  461.             for (final ForceModel forceModel : forceModels) {
  462.                 forceModel.init(initialState, target);
  463.             }
  464.         }

  465.         /** {@inheritDoc} */
  466.         @Override
  467.         public double[] computeDerivatives(final SpacecraftState state) {

  468.             currentState = state;
  469.             Arrays.fill(yDot, 0.0);
  470.             if (getOrbitType() != null) {
  471.                 // propagation uses regular orbits
  472.                 currentState.getOrbit().getJacobianWrtCartesian(getPositionAngleType(), jacobian);
  473.             }

  474.             // compute the contributions of all perturbing forces,
  475.             // using the Kepler contribution at the end since
  476.             // NewtonianAttraction is always the last instance in the list
  477.             for (final ForceModel forceModel : forceModels) {
  478.                 forceModel.addContribution(state, this);
  479.             }

  480.             if (getOrbitType() == null) {
  481.                 // position derivative is velocity, and was not added above in the force models
  482.                 // (it is added when orbit type is non-null because NewtonianAttraction considers it)
  483.                 final Vector3D velocity = currentState.getPVCoordinates().getVelocity();
  484.                 yDot[0] += velocity.getX();
  485.                 yDot[1] += velocity.getY();
  486.                 yDot[2] += velocity.getZ();
  487.             }


  488.             return yDot.clone();

  489.         }

  490.         /** {@inheritDoc} */
  491.         @Override
  492.         public void addKeplerContribution(final double mu) {
  493.             if (getOrbitType() == null) {

  494.                 // if mu is neither 0 nor NaN, we want to include Newtonian acceleration
  495.                 if (mu > 0) {
  496.                     // velocity derivative is Newtonian acceleration
  497.                     final Vector3D position = currentState.getPVCoordinates().getPosition();
  498.                     final double r2         = position.getNormSq();
  499.                     final double coeff      = -mu / (r2 * FastMath.sqrt(r2));
  500.                     yDot[3] += coeff * position.getX();
  501.                     yDot[4] += coeff * position.getY();
  502.                     yDot[5] += coeff * position.getZ();
  503.                 }

  504.             } else {
  505.                 // propagation uses regular orbits
  506.                 currentState.getOrbit().addKeplerContribution(getPositionAngleType(), mu, yDot);
  507.             }
  508.         }

  509.         /** {@inheritDoc} */
  510.         public void addNonKeplerianAcceleration(final Vector3D gamma) {
  511.             for (int i = 0; i < 6; ++i) {
  512.                 final double[] jRow = jacobian[i];
  513.                 yDot[i] += jRow[3] * gamma.getX() + jRow[4] * gamma.getY() + jRow[5] * gamma.getZ();
  514.             }
  515.         }

  516.         /** {@inheritDoc} */
  517.         @Override
  518.         public void addMassDerivative(final double q) {
  519.             if (q > 0) {
  520.                 throw new OrekitIllegalArgumentException(OrekitMessages.POSITIVE_FLOW_RATE, q);
  521.             }
  522.             yDot[6] += q;
  523.         }

  524.     }

  525.     /** Estimate tolerance vectors for integrators when propagating in absolute position-velocity-acceleration.
  526.      * @param dP user specified position error
  527.      * @param absPva reference absolute position-velocity-acceleration
  528.      * @return a two rows array, row 0 being the absolute tolerance error and row 1
  529.      * being the relative tolerance error
  530.      * @see NumericalPropagator#tolerances(double, Orbit, OrbitType)
  531.      */
  532.     public static double[][] tolerances(final double dP, final AbsolutePVCoordinates absPva) {

  533.         final double relative = dP / absPva.getPosition().getNorm();
  534.         final double dV = relative * absPva.getVelocity().getNorm();

  535.         final double[] absTol = new double[7];
  536.         final double[] relTol = new double[7];

  537.         absTol[0] = dP;
  538.         absTol[1] = dP;
  539.         absTol[2] = dP;
  540.         absTol[3] = dV;
  541.         absTol[4] = dV;
  542.         absTol[5] = dV;

  543.         // we set the mass tolerance arbitrarily to 1.0e-6 kg, as mass evolves linearly
  544.         // with trust, this often has no influence at all on propagation
  545.         absTol[6] = 1.0e-6;

  546.         Arrays.fill(relTol, relative);

  547.         return new double[][] {
  548.             absTol, relTol
  549.         };

  550.     }

  551.     /** Estimate tolerance vectors for integrators when propagating in orbits.
  552.      * <p>
  553.      * The errors are estimated from partial derivatives properties of orbits,
  554.      * starting from a scalar position error specified by the user.
  555.      * Considering the energy conservation equation V = sqrt(mu (2/r - 1/a)),
  556.      * we get at constant energy (i.e. on a Keplerian trajectory):
  557.      * <pre>
  558.      * V² r |dV| = mu |dr|
  559.      * </pre>
  560.      * <p> So we deduce a scalar velocity error consistent with the position error.
  561.      * From here, we apply orbits Jacobians matrices to get consistent errors
  562.      * on orbital parameters.
  563.      *
  564.      * <p>
  565.      * The tolerances are only <em>orders of magnitude</em>, and integrator tolerances
  566.      * are only local estimates, not global ones. So some care must be taken when using
  567.      * these tolerances. Setting 1mm as a position error does NOT mean the tolerances
  568.      * will guarantee a 1mm error position after several orbits integration.
  569.      * </p>
  570.      * @param dP user specified position error
  571.      * @param orbit reference orbit
  572.      * @param type propagation type for the meaning of the tolerance vectors elements
  573.      * (it may be different from {@code orbit.getType()})
  574.      * @return a two rows array, row 0 being the absolute tolerance error and row 1
  575.      * being the relative tolerance error
  576.      */
  577.     public static double[][] tolerances(final double dP, final Orbit orbit, final OrbitType type) {

  578.         // estimate the scalar velocity error
  579.         final PVCoordinates pv = orbit.getPVCoordinates();
  580.         final double r2 = pv.getPosition().getNormSq();
  581.         final double v  = pv.getVelocity().getNorm();
  582.         final double dV = orbit.getMu() * dP / (v * r2);

  583.         return tolerances(dP, dV, orbit, type);

  584.     }

  585.     /** Estimate tolerance vectors for integrators when propagating in orbits.
  586.      * <p>
  587.      * The errors are estimated from partial derivatives properties of orbits,
  588.      * starting from scalar position and velocity errors specified by the user.
  589.      * <p>
  590.      * The tolerances are only <em>orders of magnitude</em>, and integrator tolerances
  591.      * are only local estimates, not global ones. So some care must be taken when using
  592.      * these tolerances. Setting 1mm as a position error does NOT mean the tolerances
  593.      * will guarantee a 1mm error position after several orbits integration.
  594.      * </p>
  595.      * @param dP user specified position error
  596.      * @param dV user specified velocity error
  597.      * @param orbit reference orbit
  598.      * @param type propagation type for the meaning of the tolerance vectors elements
  599.      * (it may be different from {@code orbit.getType()})
  600.      * @return a two rows array, row 0 being the absolute tolerance error and row 1
  601.      * being the relative tolerance error
  602.      * @since 10.3
  603.      */
  604.     public static double[][] tolerances(final double dP, final double dV,
  605.                                         final Orbit orbit, final OrbitType type) {

  606.         final double[] absTol = new double[7];
  607.         final double[] relTol = new double[7];

  608.         // we set the mass tolerance arbitrarily to 1.0e-6 kg, as mass evolves linearly
  609.         // with trust, this often has no influence at all on propagation
  610.         absTol[6] = 1.0e-6;

  611.         if (type == OrbitType.CARTESIAN) {
  612.             absTol[0] = dP;
  613.             absTol[1] = dP;
  614.             absTol[2] = dP;
  615.             absTol[3] = dV;
  616.             absTol[4] = dV;
  617.             absTol[5] = dV;
  618.         } else {

  619.             // convert the orbit to the desired type
  620.             final double[][] jacobian = new double[6][6];
  621.             final Orbit converted = type.convertType(orbit);
  622.             converted.getJacobianWrtCartesian(PositionAngle.TRUE, jacobian);

  623.             for (int i = 0; i < 6; ++i) {
  624.                 final double[] row = jacobian[i];
  625.                 absTol[i] = FastMath.abs(row[0]) * dP +
  626.                             FastMath.abs(row[1]) * dP +
  627.                             FastMath.abs(row[2]) * dP +
  628.                             FastMath.abs(row[3]) * dV +
  629.                             FastMath.abs(row[4]) * dV +
  630.                             FastMath.abs(row[5]) * dV;
  631.                 if (Double.isNaN(absTol[i])) {
  632.                     throw new OrekitException(OrekitMessages.SINGULAR_JACOBIAN_FOR_ORBIT_TYPE, type);
  633.                 }
  634.             }

  635.         }

  636.         Arrays.fill(relTol, dP / FastMath.sqrt(orbit.getPVCoordinates().getPosition().getNormSq()));

  637.         return new double[][] {
  638.             absTol, relTol
  639.         };

  640.     }

  641.     /** {@inheritDoc} */
  642.     @Override
  643.     protected void beforeIntegration(final SpacecraftState initialState, final AbsoluteDate tEnd) {

  644.         if (!getFrame().isPseudoInertial()) {

  645.             // inspect all force models to find InertialForces
  646.             for (ForceModel force : forceModels) {
  647.                 if (force instanceof InertialForces) {
  648.                     return;
  649.                 }
  650.             }

  651.             // throw exception if no inertial forces found
  652.             throw new OrekitException(OrekitMessages.INERTIAL_FORCE_MODEL_MISSING, getFrame().getName());

  653.         }

  654.     }

  655. }