FieldAuxiliaryElements.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.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.  *  @author Bryan Cazabonne
  30.  * @param <T> type of the field elements
  31.  */
  32. public class FieldAuxiliaryElements<T extends CalculusFieldElement<T>> {

  33.     /** Orbit. */
  34.     private final FieldOrbit<T> orbit;

  35.     /** Orbit date. */
  36.     private final FieldAbsoluteDate<T> date;

  37.     /** Orbit frame. */
  38.     private final Frame frame;

  39.     /** Eccentricity. */
  40.     private final T ecc;

  41.     /** Keplerian mean motion. */
  42.     private final T n;

  43.     /** Keplerian period. */
  44.     private final T period;

  45.     /** Semi-major axis. */
  46.     private final T sma;

  47.     /** x component of eccentricity vector. */
  48.     private final T k;

  49.     /** y component of eccentricity vector. */
  50.     private final T h;

  51.     /** x component of inclination vector. */
  52.     private final T q;

  53.     /** y component of inclination vector. */
  54.     private final T p;

  55.     /** Mean longitude. */
  56.     private final T lm;

  57.     /** True longitude. */
  58.     private final T lv;

  59.     /** Eccentric longitude. */
  60.     private final T le;

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

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

  77.     /** C = 1 + p² + q². */
  78.     private final T C;

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

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

  83.     /** Equinoctial frame w vector. */
  84.     private final FieldVector3D<T> w;

  85.     /** Direction cosine α. */
  86.     private final T alpha;

  87.     /** Direction cosine β. */
  88.     private final T beta;

  89.     /** Direction cosine γ. */
  90.     private final T gamma;

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

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

  97.         // Orbit
  98.         this.orbit = orbit;

  99.         // Date of the orbit
  100.         date = orbit.getDate();

  101.         // Orbit definition frame
  102.         frame = orbit.getFrame();

  103.         // Eccentricity
  104.         ecc = orbit.getE();

  105.         // Keplerian mean motion
  106.         n = orbit.getKeplerianMeanMotion();

  107.         // Keplerian period
  108.         period = orbit.getKeplerianPeriod();

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

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

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

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

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

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

  140.     /** Get the orbit.
  141.      * @return the orbit
  142.      */
  143.     public FieldOrbit<T> getOrbit() {
  144.         return orbit;
  145.     }

  146.     /** Get the date of the orbit.
  147.      * @return the date
  148.      */
  149.     public FieldAbsoluteDate<T> getDate() {
  150.         return date;
  151.     }

  152.     /** Get the definition frame of the orbit.
  153.      * @return the definition frame
  154.      */
  155.     public Frame getFrame() {
  156.         return frame;
  157.     }

  158.     /** Get the eccentricity.
  159.      * @return ecc
  160.      */
  161.     public T getEcc() {
  162.         return ecc;
  163.     }

  164.     /** Get the Keplerian mean motion.
  165.      * @return n
  166.      */
  167.     public T getMeanMotion() {
  168.         return n;
  169.     }

  170.     /** Get the Keplerian period.
  171.      * @return period
  172.      */
  173.     public T getKeplerianPeriod() {
  174.         return period;
  175.     }

  176.     /** Get the semi-major axis.
  177.      * @return the semi-major axis a
  178.      */
  179.     public T getSma() {
  180.         return sma;
  181.     }

  182.     /** Get the x component of eccentricity vector.
  183.      * <p>
  184.      * This element called k in DSST corresponds to ex
  185.      * for the {@link org.orekit.orbits.EquinoctialOrbit}
  186.      * </p>
  187.      * @return k
  188.      */
  189.     public T getK() {
  190.         return k;
  191.     }

  192.     /** Get the y component of eccentricity vector.
  193.      * <p>
  194.      * This element called h in DSST corresponds to ey
  195.      * for the {@link org.orekit.orbits.EquinoctialOrbit}
  196.      * </p>
  197.      * @return h
  198.      */
  199.     public T getH() {
  200.         return h;
  201.     }

  202.     /** Get the x component of inclination vector.
  203.      * <p>
  204.      * This element called q in DSST corresponds to hx
  205.      * for the {@link org.orekit.orbits.EquinoctialOrbit}
  206.      * </p>
  207.      * @return q
  208.      */
  209.     public T getQ() {
  210.         return q;
  211.     }

  212.     /** Get the y component of inclination vector.
  213.      * <p>
  214.      * This element called p in DSST corresponds to hy
  215.      * for the {@link org.orekit.orbits.EquinoctialOrbit}
  216.      * </p>
  217.      * @return p
  218.      */
  219.     public T getP() {
  220.         return p;
  221.     }

  222.     /** Get the mean longitude.
  223.      * @return lm
  224.      */
  225.     public T getLM() {
  226.         return lm;
  227.     }

  228.     /** Get the true longitude.
  229.      * @return lv
  230.      */
  231.     public T getLv() {
  232.         return lv;
  233.     }

  234.     /** Get the eccentric longitude.
  235.      * @return le
  236.      */
  237.     public T getLe() {
  238.         return le;
  239.     }

  240.     /** Get the retrograde factor.
  241.      * @return the retrograde factor I
  242.      */
  243.     public int getRetrogradeFactor() {
  244.         return I;
  245.     }

  246.     /** Get B = sqrt(1 - e²).
  247.      * @return B
  248.      */
  249.     public T getB() {
  250.         return B;
  251.     }

  252.     /** Get C = 1 + p² + q².
  253.      * @return C
  254.      */
  255.     public T getC() {
  256.         return C;
  257.     }

  258.     /** Get equinoctial frame vector f.
  259.      * @return f vector
  260.      */
  261.     public FieldVector3D<T> getVectorF() {
  262.         return f;
  263.     }

  264.     /** Get equinoctial frame vector g.
  265.      * @return g vector
  266.      */
  267.     public FieldVector3D<T> getVectorG() {
  268.         return g;
  269.     }

  270.     /** Get equinoctial frame vector w.
  271.      * @return w vector
  272.      */
  273.     public FieldVector3D<T> getVectorW() {
  274.         return w;
  275.     }

  276.     /** Get direction cosine α for central body.
  277.      * @return α
  278.      */
  279.     public T getAlpha() {
  280.         return alpha;
  281.     }

  282.     /** Get direction cosine β for central body.
  283.      * @return β
  284.      */
  285.     public T getBeta() {
  286.         return beta;
  287.     }

  288.     /** Get direction cosine γ for central body.
  289.      * @return γ
  290.      */
  291.     public T getGamma() {
  292.         return gamma;
  293.     }

  294.     /** Transforms the FieldAuxiliaryElements instance into an AuxiliaryElements instance.
  295.      * @return the AuxiliaryElements instance
  296.      * @since 11.3.3
  297.      */
  298.     public AuxiliaryElements toAuxiliaryElements() {
  299.         return new AuxiliaryElements(orbit.toOrbit(), getRetrogradeFactor());
  300.     }
  301. }