AuxiliaryElements.java

  1. /* Copyright 2002-2013 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.propagation.semianalytical.dsst.utilities;

  18. import org.apache.commons.math3.geometry.euclidean.threed.Vector3D;
  19. import org.apache.commons.math3.util.FastMath;
  20. import org.apache.commons.math3.util.MathUtils;
  21. import org.orekit.frames.Frame;
  22. import org.orekit.orbits.Orbit;
  23. import org.orekit.time.AbsoluteDate;


  24. /** Container class for common parameters used by all DSST forces.
  25.  *  <p>
  26.  *  Most of them are defined in Danielson paper at § 2.1.
  27.  *  </p>
  28.  *  @author Pascal Parraud
  29.  */
  30. public class AuxiliaryElements {

  31.     /** Orbit date. */
  32.     private final AbsoluteDate date;

  33.     /** Orbit frame. */
  34.     private final Frame frame;

  35.     /** Central body attraction coefficient. */
  36.     private final double mu;

  37.     /** Eccentricity. */
  38.     private final double ecc;

  39.     /** Keplerian mean motion. */
  40.     private final double n;

  41.     /** Keplerian period. */
  42.     private final double period;

  43.     /** Semi-major axis. */
  44.     private final double sma;

  45.     /** x component of eccentricity vector. */
  46.     private final double k;

  47.     /** y component of eccentricity vector. */
  48.     private final double h;

  49.     /** x component of inclination vector. */
  50.     private final double q;

  51.     /** y component of inclination vector. */
  52.     private final double p;

  53.     /** Mean longitude. */
  54.     private final double lm;

  55.     /** Retrograde factor I.
  56.      *  <p>
  57.      *  DSST model needs equinoctial orbit as internal representation.
  58.      *  Classical equinoctial elements have discontinuities when inclination
  59.      *  is close to zero. In this representation, I = +1. <br>
  60.      *  To avoid this discontinuity, another representation exists and equinoctial
  61.      *  elements can be expressed in a different way, called "retrograde" orbit.
  62.      *  This implies I = -1. <br>
  63.      *  As Orekit doesn't implement the retrograde orbit, I is always set to +1.
  64.      *  But for the sake of consistency with the theory, the retrograde factor
  65.      *  has been kept in the formulas.
  66.      *  </p>
  67.      */
  68.     private final int    I;

  69.     /** A = sqrt(&mu; * a). */
  70.     private final double A;

  71.     /** B = sqrt(1 - h<sup>2</sup> - k<sup>2</sup>). */
  72.     private final double B;

  73.     /** C = 1 + p<sup>2</sup> + q<sup>2</sup>. */
  74.     private final double C;

  75.     /** Equinoctial frame f vector. */
  76.     private final Vector3D f;

  77.     /** Equinoctial frame g vector. */
  78.     private final Vector3D g;

  79.     /** Equinoctial frame w vector. */
  80.     private final Vector3D w;

  81.     /** Direction cosine &alpha;. */
  82.     private final double alpha;

  83.     /** Direction cosine &beta;. */
  84.     private final double beta;

  85.     /** Direction cosine &gamma;. */
  86.     private final double gamma;

  87.     /** Simple constructor.
  88.      * @param orbit related mean orbit for auxiliary elements
  89.      * @param retrogradeFactor retrograde factor I [Eq. 2.1.2-(2)]
  90.      */
  91.     public AuxiliaryElements(final Orbit orbit, final int retrogradeFactor) {
  92.         // Date of the orbit
  93.         date = orbit.getDate();

  94.         // Orbit definition frame
  95.         frame = orbit.getFrame();

  96.         // Central body attraction coefficient
  97.         mu = orbit.getMu();

  98.         // Eccentricity
  99.         ecc = orbit.getE();

  100.         // Keplerian mean motion
  101.         n = orbit.getKeplerianMeanMotion();

  102.         // Keplerian period
  103.         period = orbit.getKeplerianPeriod();

  104.         // Equinoctial elements [Eq. 2.1.2-(1)]
  105.         sma = orbit.getA();
  106.         k   = orbit.getEquinoctialEx();
  107.         h   = orbit.getEquinoctialEy();
  108.         q   = orbit.getHx();
  109.         p   = orbit.getHy();
  110.         lm  = MathUtils.normalizeAngle(orbit.getLM(), FastMath.PI);

  111.         // Retrograde factor [Eq. 2.1.2-(2)]
  112.         I = retrogradeFactor;

  113.         final double k2 = k * k;
  114.         final double h2 = h * h;
  115.         final double q2 = q * q;
  116.         final double p2 = p * p;

  117.         // A, B, C parameters [Eq. 2.1.6-(1)]
  118.         A = FastMath.sqrt(mu * sma);
  119.         B = FastMath.sqrt(1 - k2 - h2);
  120.         C = 1 + q2 + p2;

  121.         // Equinoctial reference frame [Eq. 2.1.4-(1)]
  122.         final double ooC = 1. / C;
  123.         final double px2 = 2. * p;
  124.         final double qx2 = 2. * q;
  125.         final double pq2 = px2 * q;
  126.         f = new Vector3D(ooC, new Vector3D(1. - p2 + q2, pq2, -px2 * I));
  127.         g = new Vector3D(ooC, new Vector3D(pq2 * I, (1. + p2 - q2) * I, qx2));
  128.         w = new Vector3D(ooC, new Vector3D(px2, -qx2, (1. - p2 - q2) * I));

  129.         // Direction cosines for central body [Eq. 2.1.9-(1)]
  130.         alpha = f.getZ();
  131.         beta  = g.getZ();
  132.         gamma = w.getZ();
  133.     }

  134.     /** Get the date of the orbit.
  135.      * @return the date
  136.      */
  137.     public AbsoluteDate getDate() {
  138.         return date;
  139.     }

  140.     /** Get the definition frame of the orbit.
  141.      * @return the definition frame
  142.      */
  143.     public Frame getFrame() {
  144.         return frame;
  145.     }

  146.     /** Get the central body attraction coefficient.
  147.      * @return &mu;
  148.      */
  149.     public double getMu() {
  150.         return mu;
  151.     }

  152.     /** Get the eccentricity.
  153.      * @return ecc
  154.      */
  155.     public double getEcc() {
  156.         return ecc;
  157.     }

  158.     /** Get the Keplerian mean motion.
  159.      * @return n
  160.      */
  161.     public double getMeanMotion() {
  162.         return n;
  163.     }

  164.     /** Get the Keplerian period.
  165.      * @return period
  166.      */
  167.     public double getKeplerianPeriod() {
  168.         return period;
  169.     }

  170.     /** Get the semi-major axis.
  171.      * @return the semi-major axis a
  172.      */
  173.     public double getSma() {
  174.         return sma;
  175.     }

  176.     /** Get the x component of eccentricity vector.
  177.      * <p>
  178.      * This element called k in DSST corresponds to ex
  179.      * for the {@link org.orekit.orbits.EquinoctialOrbit}
  180.      * </p>
  181.      * @return k
  182.      */
  183.     public double getK() {
  184.         return k;
  185.     }

  186.     /** Get the y component of eccentricity vector.
  187.      * <p>
  188.      * This element called h in DSST corresponds to ey
  189.      * for the {@link org.orekit.orbits.EquinoctialOrbit}
  190.      * </p>
  191.      * @return h
  192.      */
  193.     public double getH() {
  194.         return h;
  195.     }

  196.     /** Get the x component of inclination vector.
  197.      * <p>
  198.      * This element called q in DSST corresponds to hx
  199.      * for the {@link org.orekit.orbits.EquinoctialOrbit}
  200.      * </p>
  201.      * @return q
  202.      */
  203.     public double getQ() {
  204.         return q;
  205.     }

  206.     /** Get the y component of inclination vector.
  207.      * <p>
  208.      * This element called p in DSST corresponds to hy
  209.      * for the {@link org.orekit.orbits.EquinoctialOrbit}
  210.      * </p>
  211.      * @return p
  212.      */
  213.     public double getP() {
  214.         return p;
  215.     }

  216.     /** Get the mean longitude.
  217.      * @return lm
  218.      */
  219.     public double getLM() {
  220.         return lm;
  221.     }

  222.     /** Get the retrograde factor.
  223.      * @return the retrograde factor I
  224.      */
  225.     public int getRetrogradeFactor() {
  226.         return I;
  227.     }

  228.     /** Get A = sqrt(&mu; * a).
  229.      * @return A
  230.      */
  231.     public double getA() {
  232.         return A;
  233.     }

  234.     /** Get B = sqrt(1 - e<sup>2</sup>).
  235.      * @return B
  236.      */
  237.     public double getB() {
  238.         return B;
  239.     }

  240.     /** Get C = 1 + p<sup>2</sup> + q<sup>2</sup>.
  241.      * @return C
  242.      */
  243.     public double getC() {
  244.         return C;
  245.     }

  246.     /** Get equinoctial frame vector f.
  247.      * @return f vector
  248.      */
  249.     public Vector3D getVectorF() {
  250.         return f;
  251.     }

  252.     /** Get equinoctial frame vector g.
  253.      * @return g vector
  254.      */
  255.     public Vector3D getVectorG() {
  256.         return g;
  257.     }

  258.     /** Get equinoctial frame vector w.
  259.      * @return w vector
  260.      */
  261.     public Vector3D getVectorW() {
  262.         return w;
  263.     }

  264.     /** Get direction cosine &alpha; for central body.
  265.      * @return &alpha;
  266.      */
  267.     public double getAlpha() {
  268.         return alpha;
  269.     }

  270.     /** Get direction cosine &beta; for central body.
  271.      * @return &beta;
  272.      */
  273.     public double getBeta() {
  274.         return beta;
  275.     }

  276.     /** Get direction cosine &gamma; for central body.
  277.      * @return &gamma;
  278.      */
  279.     public double getGamma() {
  280.         return gamma;
  281.     }

  282. }