FieldAuxiliaryElements.java

  1. /* Copyright 2002-2025 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.     /** Simple constructor.
  86.      * @param orbit related mean orbit for auxiliary elements
  87.      * @param retrogradeFactor retrograde factor I [Eq. 2.1.2-(2)]
  88.      */
  89.     public FieldAuxiliaryElements(final FieldOrbit<T> orbit, final int retrogradeFactor) {

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

  91.         // Orbit
  92.         this.orbit = orbit;

  93.         // Date of the orbit
  94.         date = orbit.getDate();

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

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

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

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

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

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

  114.         final T k2 = k.multiply(k);
  115.         final T h2 = h.multiply(h);
  116.         final T q2 = q.multiply(q);
  117.         final T p2 = p.multiply(p);

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

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

  130.     /** Get the orbit.
  131.      * @return the orbit
  132.      */
  133.     public FieldOrbit<T> getOrbit() {
  134.         return orbit;
  135.     }

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

  266.     /** Transforms the FieldAuxiliaryElements instance into an AuxiliaryElements instance.
  267.      * @return the AuxiliaryElements instance
  268.      * @since 11.3.3
  269.      */
  270.     public AuxiliaryElements toAuxiliaryElements() {
  271.         return new AuxiliaryElements(orbit.toOrbit(), getRetrogradeFactor());
  272.     }

  273. }