FieldAuxiliaryElements.java

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

  18. import org.hipparchus.CalculusFieldElement;
  19. import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
  20. import org.hipparchus.util.FastMath;
  21. import org.hipparchus.util.MathUtils;
  22. import org.orekit.frames.Frame;
  23. import org.orekit.orbits.FieldOrbit;
  24. import org.orekit.time.FieldAbsoluteDate;

  25. /** Container class for common parameters used by all DSST forces.
  26.  *  <p>
  27.  *  Most of them are defined in Danielson paper at § 2.1.
  28.  *  </p>
  29.  */
  30. public class FieldAuxiliaryElements<T extends CalculusFieldElement<T>> {

  31.     /** Orbit. */
  32.     private final FieldOrbit<T> orbit;

  33.     /** Orbit date. */
  34.     private final FieldAbsoluteDate<T> date;

  35.     /** Orbit frame. */
  36.     private final Frame frame;

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

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

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

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

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

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

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

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

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

  55.     /** True longitude. */
  56.     private final T lv;

  57.     /** Eccentric longitude. */
  58.     private final T le;

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

  73.     /** B = sqrt(1 - h² - k²). */
  74.     private final T B;

  75.     /** C = 1 + p² + q². */
  76.     private final T C;

  77.     /** Equinoctial frame f vector. */
  78.     private final FieldVector3D<T> f;

  79.     /** Equinoctial frame g vector. */
  80.     private final FieldVector3D<T> g;

  81.     /** Equinoctial frame w vector. */
  82.     private final FieldVector3D<T> w;

  83.     /** Direction cosine α. */
  84.     private final T alpha;

  85.     /** Direction cosine β. */
  86.     private final T beta;

  87.     /** Direction cosine γ. */
  88.     private final T gamma;

  89.     /** Simple constructor.
  90.      * @param orbit related mean orbit for auxiliary elements
  91.      * @param retrogradeFactor retrograde factor I [Eq. 2.1.2-(2)]
  92.      */
  93.     public FieldAuxiliaryElements(final FieldOrbit<T> orbit, final int retrogradeFactor) {

  94.         final T pi = orbit.getDate().getField().getZero().getPi();

  95.         // Orbit
  96.         this.orbit = orbit;

  97.         // Date of the orbit
  98.         date = orbit.getDate();

  99.         // Orbit definition frame
  100.         frame = orbit.getFrame();

  101.         // Eccentricity
  102.         ecc = orbit.getE();

  103.         // Keplerian mean motion
  104.         n = orbit.getKeplerianMeanMotion();

  105.         // Keplerian period
  106.         period = orbit.getKeplerianPeriod();

  107.         // Equinoctial elements [Eq. 2.1.2-(1)]
  108.         sma = orbit.getA();
  109.         k   = orbit.getEquinoctialEx();
  110.         h   = orbit.getEquinoctialEy();
  111.         q   = orbit.getHx();
  112.         p   = orbit.getHy();
  113.         lm  = MathUtils.normalizeAngle(orbit.getLM(), pi);
  114.         lv  = MathUtils.normalizeAngle(orbit.getLv(), pi);
  115.         le  = MathUtils.normalizeAngle(orbit.getLE(), pi);

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

  118.         final T k2 = k.multiply(k);
  119.         final T h2 = h.multiply(h);
  120.         final T q2 = q.multiply(q);
  121.         final T p2 = p.multiply(p);

  122.         // A, B, C parameters [Eq. 2.1.6-(1)]
  123.         B = FastMath.sqrt(k2.add(h2).negate().add(1.));
  124.         C = q2.add(p2).add(1);

  125.         // Equinoctial reference frame [Eq. 2.1.4-(1)]
  126.         final T ooC = C.reciprocal();
  127.         final T px2 = p.multiply(2.);
  128.         final T qx2 = q.multiply(2.);
  129.         final T pq2 = px2.multiply(q);
  130.         f = new FieldVector3D<>(ooC, new FieldVector3D<>(p2.negate().add(1.).add(q2), pq2, px2.multiply(I).negate()));
  131.         g = new FieldVector3D<>(ooC, new FieldVector3D<>(pq2.multiply(I), (p2.add(1.).subtract(q2)).multiply(I), qx2));
  132.         w = new FieldVector3D<>(ooC, new FieldVector3D<>(px2, qx2.negate(), (p2.add(q2).negate().add(1.)).multiply(I)));

  133.         // Direction cosines for central body [Eq. 2.1.9-(1)]
  134.         alpha = (T) f.getZ();
  135.         beta  = (T) g.getZ();
  136.         gamma = (T) w.getZ();
  137.     }

  138.     /** Get the date of the orbit.
  139.      * @return the date
  140.      */
  141.     public FieldAbsoluteDate<T> getDate() {
  142.         return date;
  143.     }

  144.     /** Get the definition frame of the orbit.
  145.      * @return the definition frame
  146.      */
  147.     public Frame getFrame() {
  148.         return frame;
  149.     }

  150.     /** Get the eccentricity.
  151.      * @return ecc
  152.      */
  153.     public T getEcc() {
  154.         return ecc;
  155.     }

  156.     /** Get the Keplerian mean motion.
  157.      * @return n
  158.      */
  159.     public T getMeanMotion() {
  160.         return n;
  161.     }

  162.     /** Get the Keplerian period.
  163.      * @return period
  164.      */
  165.     public T getKeplerianPeriod() {
  166.         return period;
  167.     }

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

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

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

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

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

  214.     /** Get the mean longitude.
  215.      * @return lm
  216.      */
  217.     public T getLM() {
  218.         return lm;
  219.     }

  220.     /** Get the true longitude.
  221.      * @return lv
  222.      */
  223.     public T getLv() {
  224.         return lv;
  225.     }

  226.     /** Get the eccentric longitude.
  227.      * @return le
  228.      */
  229.     public T getLe() {
  230.         return le;
  231.     }

  232.     /** Get the retrograde factor.
  233.      * @return the retrograde factor I
  234.      */
  235.     public int getRetrogradeFactor() {
  236.         return I;
  237.     }

  238.     /** Get B = sqrt(1 - e²).
  239.      * @return B
  240.      */
  241.     public T getB() {
  242.         return B;
  243.     }

  244.     /** Get C = 1 + p² + q².
  245.      * @return C
  246.      */
  247.     public T getC() {
  248.         return C;
  249.     }

  250.     /** Get equinoctial frame vector f.
  251.      * @return f vector
  252.      */
  253.     public FieldVector3D<T> getVectorF() {
  254.         return f;
  255.     }

  256.     /** Get equinoctial frame vector g.
  257.      * @return g vector
  258.      */
  259.     public FieldVector3D<T> getVectorG() {
  260.         return g;
  261.     }

  262.     /** Get equinoctial frame vector w.
  263.      * @return w vector
  264.      */
  265.     public FieldVector3D<T> getVectorW() {
  266.         return w;
  267.     }

  268.     /** Get direction cosine α for central body.
  269.      * @return α
  270.      */
  271.     public T getAlpha() {
  272.         return alpha;
  273.     }

  274.     /** Get direction cosine β for central body.
  275.      * @return β
  276.      */
  277.     public T getBeta() {
  278.         return beta;
  279.     }

  280.     /** Get direction cosine γ for central body.
  281.      * @return γ
  282.      */
  283.     public T getGamma() {
  284.         return gamma;
  285.     }

  286.     /** Transforms the FieldAuxiliaryElements instance into an AuxiliaryElements instance.
  287.      * @return the AuxiliaryElements instance
  288.      * @since 11.3.3
  289.      */
  290.     public AuxiliaryElements toAuxiliaryElements() {
  291.         return new AuxiliaryElements(orbit.toOrbit(), getRetrogradeFactor());
  292.     }
  293. }