FieldNumericalPropagator.java

  1. /* Copyright 2002-2024 CS GROUP
  2.  * Licensed to CS GROUP (CS) under one or more
  3.  * contributor license agreements.  See the NOTICE file distributed with
  4.  * this work for additional information regarding copyright ownership.
  5.  * CS licenses this file to You under the Apache License, Version 2.0
  6.  * (the "License"); you may not use this file except in compliance with
  7.  * the License.  You may obtain a copy of the License at
  8.  *
  9.  *   http://www.apache.org/licenses/LICENSE-2.0
  10.  *
  11.  * Unless required by applicable law or agreed to in writing, software
  12.  * distributed under the License is distributed on an "AS IS" BASIS,
  13.  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
  14.  * See the License for the specific language governing permissions and
  15.  * limitations under the License.
  16.  */
  17. package org.orekit.propagation.numerical;

  18. import java.util.ArrayList;
  19. import java.util.Arrays;
  20. import java.util.Collections;
  21. import java.util.List;

  22. import org.hipparchus.CalculusFieldElement;
  23. import org.hipparchus.Field;
  24. import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
  25. import org.hipparchus.ode.FieldODEIntegrator;
  26. import org.hipparchus.util.MathArrays;
  27. import org.orekit.annotation.DefaultDataContext;
  28. import org.orekit.attitudes.AttitudeProvider;
  29. import org.orekit.attitudes.FieldAttitude;
  30. import org.orekit.data.DataContext;
  31. import org.orekit.errors.OrekitException;
  32. import org.orekit.errors.OrekitIllegalArgumentException;
  33. import org.orekit.errors.OrekitInternalError;
  34. import org.orekit.errors.OrekitMessages;
  35. import org.orekit.forces.ForceModel;
  36. import org.orekit.forces.gravity.NewtonianAttraction;
  37. import org.orekit.frames.Frame;
  38. import org.orekit.orbits.FieldOrbit;
  39. import org.orekit.orbits.OrbitType;
  40. import org.orekit.orbits.PositionAngleType;
  41. import org.orekit.propagation.FieldSpacecraftState;
  42. import org.orekit.propagation.PropagationType;
  43. import org.orekit.propagation.Propagator;
  44. import org.orekit.propagation.events.FieldEventDetector;
  45. import org.orekit.propagation.integration.FieldAbstractIntegratedPropagator;
  46. import org.orekit.propagation.integration.FieldStateMapper;
  47. import org.orekit.time.AbsoluteDate;
  48. import org.orekit.time.FieldAbsoluteDate;
  49. import org.orekit.utils.FieldAbsolutePVCoordinates;
  50. import org.orekit.utils.FieldPVCoordinates;
  51. import org.orekit.utils.ParameterDriver;
  52. import org.orekit.utils.ParameterObserver;
  53. import org.orekit.utils.TimeSpanMap;
  54. import org.orekit.utils.TimeStampedFieldPVCoordinates;

  55. /** This class propagates {@link org.orekit.orbits.FieldOrbit orbits} using
  56.  * numerical integration.
  57.  * <p>Numerical propagation is much more accurate than analytical propagation
  58.  * like for example {@link org.orekit.propagation.analytical.KeplerianPropagator
  59.  * Keplerian} or {@link org.orekit.propagation.analytical.EcksteinHechlerPropagator
  60.  * Eckstein-Hechler}, but requires a few more steps to set up to be used properly.
  61.  * Whereas analytical propagators are configured only thanks to their various
  62.  * constructors and can be used immediately after construction, numerical propagators
  63.  * configuration involve setting several parameters between construction time
  64.  * and propagation time.</p>
  65.  * <p>The configuration parameters that can be set are:</p>
  66.  * <ul>
  67.  *   <li>the initial spacecraft state ({@link #setInitialState(FieldSpacecraftState)})</li>
  68.  *   <li>the central attraction coefficient ({@link #setMu(CalculusFieldElement)})</li>
  69.  *   <li>the various force models ({@link #addForceModel(ForceModel)},
  70.  *   {@link #removeForceModels()})</li>
  71.  *   <li>the {@link OrbitType type} of orbital parameters to be used for propagation
  72.  *   ({@link #setOrbitType(OrbitType)}),
  73.  *   <li>the {@link PositionAngleType type} of position angle to be used in orbital parameters
  74.  *   to be used for propagation where it is relevant ({@link
  75.  *   #setPositionAngleType(PositionAngleType)}),
  76.  *   <li>whether {@link org.orekit.propagation.integration.FieldAdditionalDerivativesProvider additional derivatives providers}
  77.  *   should be propagated along with orbital state
  78.  *   ({@link #addAdditionalDerivativesProvider(org.orekit.propagation.integration.FieldAdditionalDerivativesProvider)}),
  79.  *   <li>the discrete events that should be triggered during propagation
  80.  *   ({@link #addEventDetector(FieldEventDetector)},
  81.  *   {@link #clearEventsDetectors()})</li>
  82.  *   <li>the binding logic with the rest of the application ({@link #getMultiplexer()})</li>
  83.  * </ul>
  84.  * <p>From these configuration parameters, only the initial state is mandatory. The default
  85.  * propagation settings are in {@link OrbitType#EQUINOCTIAL equinoctial} parameters with
  86.  * {@link PositionAngleType#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.FieldEquinoctialOrbit 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.FieldKeplerianOrbit Keplerian orbit parameters} (a, e, i, ω, Ω,
  103.  *   M or E or v) in meters and radians,</li>
  104.  *   <li>the {@link org.orekit.orbits.FieldCircularOrbit 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.FieldCartesianOrbit 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.  * The last element is the mass in kilograms.
  111.  * <p>The following code snippet shows a typical setting for Low Earth Orbit propagation in
  112.  * equinoctial parameters and true longitude argument:</p>
  113.  * <pre>
  114.  * final T          zero      = field.getZero();
  115.  * final T          dP        = zero.add(0.001);
  116.  * final T          minStep   = zero.add(0.001);
  117.  * final T          maxStep   = zero.add(500);
  118.  * final T          initStep  = zero.add(60);
  119.  * final double[][] tolerance = FieldNumericalPropagator.tolerances(dP, orbit, OrbitType.EQUINOCTIAL);
  120.  * AdaptiveStepsizeFieldIntegrator&lt;T&gt; integrator = new DormandPrince853FieldIntegrator&lt;&gt;(field, minStep, maxStep, tolerance[0], tolerance[1]);
  121.  * integrator.setInitialStepSize(initStep);
  122.  * propagator = new FieldNumericalPropagator&lt;&gt;(field, 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.  * @see FieldSpacecraftState
  131.  * @see ForceModel
  132.  * @see org.orekit.propagation.sampling.FieldOrekitStepHandler
  133.  * @see org.orekit.propagation.sampling.FieldOrekitFixedStepHandler
  134.  * @see org.orekit.propagation.integration.FieldIntegratedEphemeris
  135.  * @see FieldTimeDerivativesEquations
  136.  *
  137.  * @author Mathieu Rom&eacute;ro
  138.  * @author Luc Maisonobe
  139.  * @author Guylaine Prat
  140.  * @author Fabien Maussion
  141.  * @author V&eacute;ronique Pommier-Maurussane
  142.  * @param <T> type of the field elements
  143.  */
  144. public class FieldNumericalPropagator<T extends CalculusFieldElement<T>> extends FieldAbstractIntegratedPropagator<T> {

  145.     /** Force models used during the extrapolation of the FieldOrbit<T>, without Jacobians. */
  146.     private final List<ForceModel> forceModels;

  147.     /** Field used by this class.*/
  148.     private final Field<T> field;

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

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

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

  198.     /** Set the flag to ignore or not the creation of a {@link NewtonianAttraction}.
  199.      * @param ignoreCentralAttraction if true, {@link NewtonianAttraction} is <em>not</em>
  200.      * added automatically if missing
  201.      */
  202.     public void setIgnoreCentralAttraction(final boolean ignoreCentralAttraction) {
  203.         this.ignoreCentralAttraction = ignoreCentralAttraction;
  204.     }

  205.      /** Set the central attraction coefficient μ.
  206.       * <p>
  207.       * Setting the central attraction coefficient is
  208.       * equivalent to {@link #addForceModel(ForceModel) add}
  209.       * a {@link NewtonianAttraction} force model.
  210.       * </p>
  211.      * @param mu central attraction coefficient (m³/s²)
  212.      * @see #addForceModel(ForceModel)
  213.      * @see #getAllForceModels()
  214.      */
  215.     public void setMu(final T mu) {
  216.         if (ignoreCentralAttraction) {
  217.             superSetMu(mu);
  218.         } else {
  219.             addForceModel(new NewtonianAttraction(mu.getReal()));
  220.         }
  221.     }

  222.     /** Set the central attraction coefficient μ only in upper class.
  223.      * @param mu central attraction coefficient (m³/s²)
  224.      */
  225.     private void superSetMu(final T mu) {
  226.         super.setMu(mu);
  227.     }

  228.     /** Check if Newtonian attraction force model is available.
  229.      * <p>
  230.      * Newtonian attraction is always the last force model in the list.
  231.      * </p>
  232.      * @return true if Newtonian attraction force model is available
  233.      */
  234.     private boolean hasNewtonianAttraction() {
  235.         final int last = forceModels.size() - 1;
  236.         return last >= 0 && forceModels.get(last) instanceof NewtonianAttraction;
  237.     }

  238.     /** Add a force model to the global perturbation model.
  239.      * <p>If this method is not called at all, the integrated orbit will follow
  240.      * a Keplerian evolution only.</p>
  241.      * @param model perturbing {@link ForceModel} to add
  242.      * @see #removeForceModels()
  243.      * @see #setMu(CalculusFieldElement)
  244.      */
  245.     public void addForceModel(final ForceModel model) {

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

  248.             try {
  249.                 // ensure we are notified of any mu change
  250.                 model.getParametersDrivers().get(0).addObserver(new ParameterObserver() {
  251.                     /** {@inheritDoc} */
  252.                     @Override
  253.                     public void valueChanged(final double previousValue, final ParameterDriver driver, final AbsoluteDate date) {
  254.                         // mu PDriver should have only 1 span
  255.                         superSetMu(field.getZero().add(driver.getValue(date)));
  256.                     }
  257.                     /** {@inheritDoc} */
  258.                     @Override
  259.                     public void valueSpanMapChanged(final TimeSpanMap<Double> previousValue, final ParameterDriver driver) {
  260.                         // mu PDriver should have only 1 span
  261.                         superSetMu(field.getZero().add(driver.getValue()));
  262.                     }
  263.                 });
  264.             } catch (OrekitException oe) {
  265.                 // this should never happen
  266.                 throw new OrekitInternalError(oe);
  267.             }

  268.             if (hasNewtonianAttraction()) {
  269.                 // there is already a central attraction model, replace it
  270.                 forceModels.set(forceModels.size() - 1, model);
  271.             } else {
  272.                 // there are no central attraction model yet, add it at the end of the list
  273.                 forceModels.add(model);
  274.             }
  275.         } else {
  276.             // we want to add a perturbing force model
  277.             if (hasNewtonianAttraction()) {
  278.                 // insert the new force model before Newtonian attraction,
  279.                 // which should always be the last one in the list
  280.                 forceModels.add(forceModels.size() - 1, model);
  281.             } else {
  282.                 // we only have perturbing force models up to now, just append at the end of the list
  283.                 forceModels.add(model);
  284.             }
  285.         }

  286.     }

  287.     /** Remove all perturbing force models from the global perturbation model.
  288.      * <p>Once all perturbing forces have been removed (and as long as no new force
  289.      * model is added), the integrated orbit will follow a Keplerian evolution
  290.      * only.</p>
  291.      * @see #addForceModel(ForceModel)
  292.      */
  293.     public void removeForceModels() {
  294.         forceModels.clear();
  295.     }

  296.     /** Get all the force models, perturbing forces and Newtonian attraction included.
  297.      * @return list of perturbing force models, with Newtonian attraction being the
  298.      * last one
  299.      * @see #addForceModel(ForceModel)
  300.      * @see #setMu(CalculusFieldElement)
  301.      * @since 9.1
  302.      */
  303.     public List<ForceModel> getAllForceModels() {
  304.         return Collections.unmodifiableList(forceModels);
  305.     }

  306.     /** Set propagation orbit type.
  307.      * @param orbitType orbit type to use for propagation
  308.      */
  309.     public void setOrbitType(final OrbitType orbitType) {
  310.         super.setOrbitType(orbitType);
  311.     }

  312.     /** Get propagation parameter type.
  313.      * @return orbit type used for propagation
  314.      */
  315.     public OrbitType getOrbitType() {
  316.         return superGetOrbitType();
  317.     }

  318.     /** Get propagation parameter type.
  319.      * @return orbit type used for propagation
  320.      */
  321.     private OrbitType superGetOrbitType() {
  322.         return super.getOrbitType();
  323.     }

  324.     /** Set position angle type.
  325.      * <p>
  326.      * The position parameter type is meaningful only if {@link
  327.      * #getOrbitType() propagation orbit type}
  328.      * support it. As an example, it is not meaningful for propagation
  329.      * in {@link OrbitType#CARTESIAN Cartesian} parameters.
  330.      * </p>
  331.      * @param positionAngleType angle type to use for propagation
  332.      */
  333.     public void setPositionAngleType(final PositionAngleType positionAngleType) {
  334.         super.setPositionAngleType(positionAngleType);
  335.     }

  336.     /** Get propagation parameter type.
  337.      * @return angle type to use for propagation
  338.      */
  339.     public PositionAngleType getPositionAngleType() {
  340.         return super.getPositionAngleType();
  341.     }

  342.     /** Set the initial state.
  343.      * @param initialState initial state
  344.      */
  345.     public void setInitialState(final FieldSpacecraftState<T> initialState) {
  346.         resetInitialState(initialState);
  347.     }

  348.     /** {@inheritDoc} */
  349.     public void resetInitialState(final FieldSpacecraftState<T> state) {
  350.         super.resetInitialState(state);
  351.         if (!hasNewtonianAttraction()) {
  352.             setMu(state.getMu());
  353.         }
  354.         setStartDate(state.getDate());
  355.     }

  356.     /** {@inheritDoc} */
  357.     public TimeStampedFieldPVCoordinates<T> getPVCoordinates(final FieldAbsoluteDate<T> date, final Frame frame) {
  358.         return propagate(date).getPVCoordinates(frame);
  359.     }

  360.     /** {@inheritDoc} */
  361.     protected FieldStateMapper<T> createMapper(final FieldAbsoluteDate<T> referenceDate, final T mu,
  362.                                        final OrbitType orbitType, final PositionAngleType positionAngleType,
  363.                                        final AttitudeProvider attitudeProvider, final Frame frame) {
  364.         return new FieldOsculatingMapper(referenceDate, mu, orbitType, positionAngleType, attitudeProvider, frame);
  365.     }

  366.     /** Internal mapper using directly osculating parameters. */
  367.     private class FieldOsculatingMapper extends FieldStateMapper<T> {

  368.         /** Simple constructor.
  369.          * <p>
  370.          * The position parameter type is meaningful only if {@link
  371.          * #getOrbitType() propagation orbit type}
  372.          * support it. As an example, it is not meaningful for propagation
  373.          * in {@link OrbitType#CARTESIAN Cartesian} parameters.
  374.          * </p>
  375.          * @param referenceDate reference date
  376.          * @param mu central attraction coefficient (m³/s²)
  377.          * @param orbitType orbit type to use for mapping
  378.          * @param positionAngleType angle type to use for propagation
  379.          * @param attitudeProvider attitude provider
  380.          * @param frame inertial frame
  381.          */
  382.         FieldOsculatingMapper(final FieldAbsoluteDate<T> referenceDate, final T mu,
  383.                               final OrbitType orbitType, final PositionAngleType positionAngleType,
  384.                               final AttitudeProvider attitudeProvider, final Frame frame) {
  385.             super(referenceDate, mu, orbitType, positionAngleType, attitudeProvider, frame);
  386.         }

  387.         /** {@inheritDoc} */
  388.         public FieldSpacecraftState<T> mapArrayToState(final FieldAbsoluteDate<T> date, final T[] y, final T[] yDot,
  389.                                                        final PropagationType type) {
  390.             // the parameter type is ignored for the Numerical Propagator

  391.             final T mass = y[6];
  392.             if (mass.getReal() <= 0.0) {
  393.                 throw new OrekitException(OrekitMessages.NOT_POSITIVE_SPACECRAFT_MASS, mass);
  394.             }

  395.             if (superGetOrbitType() == null) {
  396.                 // propagation uses absolute position-velocity-acceleration
  397.                 final FieldVector3D<T> p = new FieldVector3D<>(y[0],    y[1],    y[2]);
  398.                 final FieldVector3D<T> v = new FieldVector3D<>(y[3],    y[4],    y[5]);
  399.                 final FieldVector3D<T> a;
  400.                 final FieldAbsolutePVCoordinates<T> absPva;
  401.                 if (yDot == null) {
  402.                     absPva = new FieldAbsolutePVCoordinates<>(getFrame(), new TimeStampedFieldPVCoordinates<>(date, p, v, FieldVector3D.getZero(date.getField())));
  403.                 } else {
  404.                     a = new FieldVector3D<>(yDot[3], yDot[4], yDot[5]);
  405.                     absPva = new FieldAbsolutePVCoordinates<>(getFrame(), new TimeStampedFieldPVCoordinates<>(date, p, v, a));
  406.                 }

  407.                 final FieldAttitude<T> attitude = getAttitudeProvider().getAttitude(absPva, date, getFrame());
  408.                 return new FieldSpacecraftState<>(absPva, attitude, mass);
  409.             } else {
  410.                 // propagation uses regular orbits
  411.                 final FieldOrbit<T> orbit       = superGetOrbitType().mapArrayToOrbit(y, yDot, super.getPositionAngleType(), date, getMu(), getFrame());
  412.                 final FieldAttitude<T> attitude = getAttitudeProvider().getAttitude(orbit, date, getFrame());
  413.                 return new FieldSpacecraftState<>(orbit, attitude, mass);
  414.             }
  415.         }

  416.         /** {@inheritDoc} */
  417.         public void mapStateToArray(final FieldSpacecraftState<T> state, final T[] y, final T[] yDot) {
  418.             if (superGetOrbitType() == null) {
  419.                 // propagation uses absolute position-velocity-acceleration
  420.                 final FieldVector3D<T> p = state.getAbsPVA().getPosition();
  421.                 final FieldVector3D<T> v = state.getAbsPVA().getVelocity();
  422.                 y[0] = p.getX();
  423.                 y[1] = p.getY();
  424.                 y[2] = p.getZ();
  425.                 y[3] = v.getX();
  426.                 y[4] = v.getY();
  427.                 y[5] = v.getZ();
  428.                 y[6] = state.getMass();
  429.             }
  430.             else {
  431.                 superGetOrbitType().mapOrbitToArray(state.getOrbit(), super.getPositionAngleType(), y, yDot);
  432.                 y[6] = state.getMass();
  433.             }
  434.         }

  435.     }

  436.     /** {@inheritDoc} */
  437.     protected MainStateEquations<T> getMainStateEquations(final FieldODEIntegrator<T> integrator) {
  438.         return new Main(integrator);
  439.     }

  440.     /** Internal class for osculating parameters integration. */
  441.     private class Main implements MainStateEquations<T>, FieldTimeDerivativesEquations<T> {

  442.         /** Derivatives array. */
  443.         private final T[] yDot;

  444.         /** Current state. */
  445.         private FieldSpacecraftState<T> currentState;

  446.         /** Jacobian of the orbital parameters with respect to the Cartesian parameters. */
  447.         private T[][] jacobian;

  448.         /** Simple constructor.
  449.          * @param integrator numerical integrator to use for propagation.
  450.          */
  451.         Main(final FieldODEIntegrator<T> integrator) {

  452.             this.yDot     = MathArrays.buildArray(getField(),  7);
  453.             this.jacobian = MathArrays.buildArray(getField(),  6, 6);
  454.             for (final ForceModel forceModel : forceModels) {
  455.                 forceModel.getFieldEventDetectors(getField()).forEach(detector -> setUpEventDetector(integrator, detector));
  456.             }

  457.             if (superGetOrbitType() == null) {
  458.                 // propagation uses absolute position-velocity-acceleration
  459.                 // we can set Jacobian once and for all
  460.                 for (int i = 0; i < jacobian.length; ++i) {
  461.                     Arrays.fill(jacobian[i], getField().getZero());
  462.                     jacobian[i][i] = getField().getOne();
  463.                 }
  464.             }

  465.         }

  466.         /** {@inheritDoc} */
  467.         @Override
  468.         public void init(final FieldSpacecraftState<T> initialState, final FieldAbsoluteDate<T> target) {
  469.             forceModels.forEach(fm -> fm.init(initialState, target));
  470.         }

  471.         /** {@inheritDoc} */
  472.         @Override
  473.         public T[] computeDerivatives(final FieldSpacecraftState<T> state) {
  474.             final T zero = state.getA().getField().getZero();
  475.             currentState = state;
  476.             Arrays.fill(yDot, zero);
  477.             if (superGetOrbitType() != null) {
  478.                 // propagation uses regular orbits
  479.                 currentState.getOrbit().getJacobianWrtCartesian(getPositionAngleType(), jacobian);
  480.             }

  481.             // compute the contributions of all perturbing forces,
  482.             // using the Kepler contribution at the end since
  483.             // NewtonianAttraction is always the last instance in the list
  484.             for (final ForceModel forceModel : forceModels) {
  485.                 forceModel.addContribution(state, this);
  486.             }

  487.             if (superGetOrbitType() == null) {
  488.                 // position derivative is velocity, and was not added above in the force models
  489.                 // (it is added when orbit type is non-null because NewtonianAttraction considers it)
  490.                 final FieldVector3D<T> velocity = currentState.getPVCoordinates().getVelocity();
  491.                 yDot[0] = yDot[0].add(velocity.getX());
  492.                 yDot[1] = yDot[1].add(velocity.getY());
  493.                 yDot[2] = yDot[2].add(velocity.getZ());
  494.             }

  495.             return yDot.clone();

  496.         }

  497.         /** {@inheritDoc} */
  498.         @Override
  499.         public void addKeplerContribution(final T mu) {
  500.             if (superGetOrbitType() == null) {

  501.                 // if mu is neither 0 nor NaN, we want to include Newtonian acceleration
  502.                 if (mu.getReal() > 0) {
  503.                     // velocity derivative is Newtonian acceleration
  504.                     final FieldVector3D<T> position = currentState.getPosition();
  505.                     final T r2         = position.getNormSq();
  506.                     final T coeff      = r2.multiply(r2.sqrt()).reciprocal().negate().multiply(mu);
  507.                     yDot[3] = yDot[3].add(coeff.multiply(position.getX()));
  508.                     yDot[4] = yDot[4].add(coeff.multiply(position.getY()));
  509.                     yDot[5] = yDot[5].add(coeff.multiply(position.getZ()));
  510.                 }

  511.             } else {
  512.                 // propagation uses regular orbits
  513.                 currentState.getOrbit().addKeplerContribution(getPositionAngleType(), mu, yDot);
  514.             }
  515.         }

  516.         /** {@inheritDoc} */
  517.         @Override
  518.         public void addNonKeplerianAcceleration(final FieldVector3D<T> gamma) {
  519.             for (int i = 0; i < 6; ++i) {
  520.                 final T[] jRow = jacobian[i];
  521.                 yDot[i] = yDot[i].add(jRow[3].linearCombination(jRow[3], gamma.getX(),
  522.                                                                 jRow[4], gamma.getY(),
  523.                                                                 jRow[5], gamma.getZ()));
  524.             }
  525.         }

  526.         /** {@inheritDoc} */
  527.         @Override
  528.         public void addMassDerivative(final T q) {
  529.             if (q.getReal() > 0) {
  530.                 throw new OrekitIllegalArgumentException(OrekitMessages.POSITIVE_FLOW_RATE, q);
  531.             }
  532.             yDot[6] = yDot[6].add(q);
  533.         }

  534.     }

  535.     /** Estimate tolerance vectors for integrators.
  536.      * <p>
  537.      * The errors are estimated from partial derivatives properties of orbits,
  538.      * starting from a scalar position error specified by the user.
  539.      * Considering the energy conservation equation V = sqrt(mu (2/r - 1/a)),
  540.      * we get at constant energy (i.e. on a Keplerian trajectory):
  541.      * <pre>
  542.      * V r² |dV| = mu |dr|
  543.      * </pre>
  544.      * So we deduce a scalar velocity error consistent with the position error.
  545.      * From here, we apply orbits Jacobians matrices to get consistent errors
  546.      * on orbital parameters.
  547.      * <p>
  548.      * The tolerances are only <em>orders of magnitude</em>, and integrator tolerances
  549.      * are only local estimates, not global ones. So some care must be taken when using
  550.      * these tolerances. Setting 1mm as a position error does NOT mean the tolerances
  551.      * will guarantee a 1mm error position after several orbits integration.
  552.      * </p>
  553.      * @param dP user specified position error
  554.      * @param orbit reference orbit
  555.      * @param type propagation type for the meaning of the tolerance vectors elements
  556.      * (it may be different from {@code orbit.getType()})
  557.      * @return a two rows array, row 0 being the absolute tolerance error and row 1
  558.      * being the relative tolerance error
  559.      * @param <T> elements type
  560.      */
  561.     public static <T extends CalculusFieldElement<T>> double[][] tolerances(final T dP, final FieldOrbit<T> orbit, final OrbitType type) {

  562.         // estimate the scalar velocity error
  563.         final FieldPVCoordinates<T> pv = orbit.getPVCoordinates();
  564.         final T r2 = pv.getPosition().getNormSq();
  565.         final T v  = pv.getVelocity().getNorm();
  566.         final T dV = dP.multiply(orbit.getMu()).divide(v.multiply(r2));

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

  568.     }

  569.     /** Estimate tolerance vectors for integrators when propagating in orbits.
  570.      * <p>
  571.      * The errors are estimated from partial derivatives properties of orbits,
  572.      * starting from scalar position and velocity errors specified by the user.
  573.      * <p>
  574.      * The tolerances are only <em>orders of magnitude</em>, and integrator tolerances
  575.      * are only local estimates, not global ones. So some care must be taken when using
  576.      * these tolerances. Setting 1mm as a position error does NOT mean the tolerances
  577.      * will guarantee a 1mm error position after several orbits integration.
  578.      * </p>
  579.      * @param <T> elements type
  580.      * @param dP user specified position error
  581.      * @param dV user specified velocity error
  582.      * @param orbit reference orbit
  583.      * @param type propagation type for the meaning of the tolerance vectors elements
  584.      * (it may be different from {@code orbit.getType()})
  585.      * @return a two rows array, row 0 being the absolute tolerance error and row 1
  586.      * being the relative tolerance error
  587.      * @since 10.3
  588.      */
  589.     public static <T extends CalculusFieldElement<T>> double[][] tolerances(final T dP, final T dV,
  590.                                                                         final FieldOrbit<T> orbit, final OrbitType type) {

  591.         final double[] absTol = new double[7];
  592.         final double[] relTol = new double[7];

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

  596.         if (type == OrbitType.CARTESIAN) {
  597.             absTol[0] = dP.getReal();
  598.             absTol[1] = dP.getReal();
  599.             absTol[2] = dP.getReal();
  600.             absTol[3] = dV.getReal();
  601.             absTol[4] = dV.getReal();
  602.             absTol[5] = dV.getReal();
  603.         } else {

  604.             // convert the orbit to the desired type
  605.             final T[][] jacobian = MathArrays.buildArray(dP.getField(), 6, 6);
  606.             final FieldOrbit<T> converted = type.convertType(orbit);
  607.             converted.getJacobianWrtCartesian(PositionAngleType.TRUE, jacobian);

  608.             for (int i = 0; i < 6; ++i) {
  609.                 final  T[] row = jacobian[i];
  610.                 absTol[i] =     row[0].abs().multiply(dP).
  611.                             add(row[1].abs().multiply(dP)).
  612.                             add(row[2].abs().multiply(dP)).
  613.                             add(row[3].abs().multiply(dV)).
  614.                             add(row[4].abs().multiply(dV)).
  615.                             add(row[5].abs().multiply(dV)).
  616.                             getReal();
  617.                 if (Double.isNaN(absTol[i])) {
  618.                     throw new OrekitException(OrekitMessages.SINGULAR_JACOBIAN_FOR_ORBIT_TYPE, type);
  619.                 }
  620.             }

  621.         }

  622.         Arrays.fill(relTol, dP.divide(orbit.getPosition().getNormSq().sqrt()).getReal());

  623.         return new double[][] { absTol, relTol };

  624.     }

  625. }