PosVelChebyshev.java

  1. /* Copyright 2002-2018 CS Systèmes d'Information
  2.  * Licensed to CS Systèmes d'Information (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.bodies;

  18. import java.io.Serializable;

  19. import org.hipparchus.RealFieldElement;
  20. import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
  21. import org.hipparchus.geometry.euclidean.threed.Vector3D;
  22. import org.orekit.time.AbsoluteDate;
  23. import org.orekit.time.FieldAbsoluteDate;
  24. import org.orekit.time.TimeScale;
  25. import org.orekit.time.TimeStamped;
  26. import org.orekit.utils.FieldPVCoordinates;
  27. import org.orekit.utils.PVCoordinates;


  28. /** Position-Velocity model based on Chebyshev polynomials.
  29.  * <p>This class represent the most basic element of the piecewise ephemerides
  30.  * for solar system bodies like JPL DE 405 ephemerides.</p>
  31.  * @see JPLEphemeridesLoader
  32.  * @author Luc Maisonobe
  33.  */
  34. class PosVelChebyshev implements TimeStamped, Serializable {

  35.     /** Serializable UID. */
  36.     private static final long serialVersionUID = 20151023L;

  37.     /** Time scale in which the ephemeris is defined. */
  38.     private final TimeScale timeScale;

  39.     /** Start of the validity range of the instance. */
  40.     private final AbsoluteDate start;

  41.     /** Duration of validity range of the instance. */
  42.     private final double duration;

  43.     /** Chebyshev polynomials coefficients for the X component. */
  44.     private final double[] xCoeffs;

  45.     /** Chebyshev polynomials coefficients for the Y component. */
  46.     private final double[] yCoeffs;

  47.     /** Chebyshev polynomials coefficients for the Z component. */
  48.     private final double[] zCoeffs;

  49.     /** Simple constructor.
  50.      * @param start start of the validity range of the instance
  51.      * @param timeScale time scale in which the ephemeris is defined
  52.      * @param duration duration of the validity range of the instance
  53.      * @param xCoeffs Chebyshev polynomials coefficients for the X component
  54.      * (a reference to the array will be stored in the instance)
  55.      * @param yCoeffs Chebyshev polynomials coefficients for the Y component
  56.      * (a reference to the array will be stored in the instance)
  57.      * @param zCoeffs Chebyshev polynomials coefficients for the Z component
  58.      * (a reference to the array will be stored in the instance)
  59.      */
  60.     PosVelChebyshev(final AbsoluteDate start, final TimeScale timeScale, final double duration,
  61.                     final double[] xCoeffs, final double[] yCoeffs, final double[] zCoeffs) {
  62.         this.start     = start;
  63.         this.timeScale = timeScale;
  64.         this.duration  = duration;
  65.         this.xCoeffs   = xCoeffs;
  66.         this.yCoeffs   = yCoeffs;
  67.         this.zCoeffs   = zCoeffs;
  68.     }

  69.     /** {@inheritDoc} */
  70.     public AbsoluteDate getDate() {
  71.         return start;
  72.     }

  73.     /** Check if a date is in validity range.
  74.      * @param date date to check
  75.      * @return true if date is in validity range
  76.      */
  77.     public boolean inRange(final AbsoluteDate date) {
  78.         final double dt = date.offsetFrom(start, timeScale);
  79.         return (dt >= -0.001) && (dt <= duration + 0.001);
  80.     }

  81.     /** Get the position-velocity-acceleration at a specified date.
  82.      * @param date date at which position-velocity-acceleration is requested
  83.      * @return position-velocity-acceleration at specified date
  84.      */
  85.     public PVCoordinates getPositionVelocityAcceleration(final AbsoluteDate date) {

  86.         // normalize date
  87.         final double t = (2 * date.offsetFrom(start, timeScale) - duration) / duration;
  88.         final double twoT = 2 * t;

  89.         // initialize Chebyshev polynomials recursion
  90.         double pKm1 = 1;
  91.         double pK   = t;
  92.         double xP   = xCoeffs[0];
  93.         double yP   = yCoeffs[0];
  94.         double zP   = zCoeffs[0];

  95.         // initialize Chebyshev polynomials derivatives recursion
  96.         double qKm1 = 0;
  97.         double qK   = 1;
  98.         double xV   = 0;
  99.         double yV   = 0;
  100.         double zV   = 0;

  101.         // initialize Chebyshev polynomials second derivatives recursion
  102.         double rKm1 = 0;
  103.         double rK   = 0;
  104.         double xA   = 0;
  105.         double yA   = 0;
  106.         double zA   = 0;

  107.         // combine polynomials by applying coefficients
  108.         for (int k = 1; k < xCoeffs.length; ++k) {

  109.             // consider last computed polynomials on position
  110.             xP += xCoeffs[k] * pK;
  111.             yP += yCoeffs[k] * pK;
  112.             zP += zCoeffs[k] * pK;

  113.             // consider last computed polynomials on velocity
  114.             xV += xCoeffs[k] * qK;
  115.             yV += yCoeffs[k] * qK;
  116.             zV += zCoeffs[k] * qK;

  117.             // consider last computed polynomials on acceleration
  118.             xA += xCoeffs[k] * rK;
  119.             yA += yCoeffs[k] * rK;
  120.             zA += zCoeffs[k] * rK;

  121.             // compute next Chebyshev polynomial value
  122.             final double pKm2 = pKm1;
  123.             pKm1 = pK;
  124.             pK   = twoT * pKm1 - pKm2;

  125.             // compute next Chebyshev polynomial derivative
  126.             final double qKm2 = qKm1;
  127.             qKm1 = qK;
  128.             qK   = twoT * qKm1 + 2 * pKm1 - qKm2;

  129.             // compute next Chebyshev polynomial second derivative
  130.             final double rKm2 = rKm1;
  131.             rKm1 = rK;
  132.             rK   = twoT * rKm1 + 4 * qKm1 - rKm2;

  133.         }

  134.         final double vScale = 2 / duration;
  135.         final double aScale = vScale * vScale;
  136.         return new PVCoordinates(new Vector3D(xP, yP, zP),
  137.                                  new Vector3D(xV * vScale, yV * vScale, zV * vScale),
  138.                                  new Vector3D(xA * aScale, yA * aScale, zA * aScale));

  139.     }

  140.     /** Get the position-velocity-acceleration at a specified date.
  141.      * @param date date at which position-velocity-acceleration is requested
  142.      * @param <T> type fo the field elements
  143.      * @return position-velocity-acceleration at specified date
  144.      */
  145.     public <T extends RealFieldElement<T>> FieldPVCoordinates<T> getPositionVelocityAcceleration(final FieldAbsoluteDate<T> date) {

  146.         final T zero = date.getField().getZero();
  147.         final T one  = date.getField().getOne();

  148.         // normalize date
  149.         final T t = date.offsetFrom(new FieldAbsoluteDate<>(date.getField(), start), timeScale).multiply(2).subtract(duration).divide(duration);
  150.         final T twoT = t.add(t);

  151.         // initialize Chebyshev polynomials recursion
  152.         T pKm1 = one;
  153.         T pK   = t;
  154.         T xP   = zero.add(xCoeffs[0]);
  155.         T yP   = zero.add(yCoeffs[0]);
  156.         T zP   = zero.add(zCoeffs[0]);

  157.         // initialize Chebyshev polynomials derivatives recursion
  158.         T qKm1 = zero;
  159.         T qK   = one;
  160.         T xV   = zero;
  161.         T yV   = zero;
  162.         T zV   = zero;

  163.         // initialize Chebyshev polynomials second derivatives recursion
  164.         T rKm1 = zero;
  165.         T rK   = zero;
  166.         T xA   = zero;
  167.         T yA   = zero;
  168.         T zA   = zero;

  169.         // combine polynomials by applying coefficients
  170.         for (int k = 1; k < xCoeffs.length; ++k) {

  171.             // consider last computed polynomials on position
  172.             xP = xP.add(pK.multiply(xCoeffs[k]));
  173.             yP = yP.add(pK.multiply(yCoeffs[k]));
  174.             zP = zP.add(pK.multiply(zCoeffs[k]));

  175.             // consider last computed polynomials on velocity
  176.             xV = xV.add(qK.multiply(xCoeffs[k]));
  177.             yV = yV.add(qK.multiply(yCoeffs[k]));
  178.             zV = zV.add(qK.multiply(zCoeffs[k]));

  179.             // consider last computed polynomials on acceleration
  180.             xA = xA.add(rK.multiply(xCoeffs[k]));
  181.             yA = yA.add(rK.multiply(yCoeffs[k]));
  182.             zA = zA.add(rK.multiply(zCoeffs[k]));

  183.             // compute next Chebyshev polynomial value
  184.             final T pKm2 = pKm1;
  185.             pKm1 = pK;
  186.             pK   = twoT.multiply(pKm1).subtract(pKm2);

  187.             // compute next Chebyshev polynomial derivative
  188.             final T qKm2 = qKm1;
  189.             qKm1 = qK;
  190.             qK   = twoT.multiply(qKm1).add(pKm1.multiply(2)).subtract(qKm2);

  191.             // compute next Chebyshev polynomial second derivative
  192.             final T rKm2 = rKm1;
  193.             rKm1 = rK;
  194.             rK   = twoT.multiply(rKm1).add(qKm1.multiply(4)).subtract(rKm2);

  195.         }

  196.         final double vScale = 2 / duration;
  197.         final double aScale = vScale * vScale;
  198.         return new FieldPVCoordinates<>(new FieldVector3D<>(xP, yP, zP),
  199.                                         new FieldVector3D<>(xV.multiply(vScale), yV.multiply(vScale), zV.multiply(vScale)),
  200.                                         new FieldVector3D<>(xA.multiply(aScale), yA.multiply(aScale), zA.multiply(aScale)));

  201.     }

  202. }