FieldDSSTGravityContext.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.forces;

  18. import org.hipparchus.CalculusFieldElement;
  19. import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
  20. import org.hipparchus.geometry.euclidean.threed.Vector3D;
  21. import org.hipparchus.util.FastMath;
  22. import org.orekit.forces.gravity.potential.UnnormalizedSphericalHarmonicsProvider;
  23. import org.orekit.frames.FieldStaticTransform;
  24. import org.orekit.frames.Frame;
  25. import org.orekit.propagation.semianalytical.dsst.utilities.FieldAuxiliaryElements;

  26. /**
  27.  * This class is a container for the common parameters used in {@link DSSTTesseral} and {@link DSSTZonal}.
  28.  * <p>
  29.  * It performs parameters initialization at each integration step for the Tesseral and Zonal contribution
  30.  * to the central body gravitational perturbation.
  31.  * </p>
  32.  * @author Bryan Cazabonne
  33.  * @author Maxime Journot
  34.  * @since 12.2
  35.  */
  36. public class FieldDSSTGravityContext<T extends CalculusFieldElement<T>> extends FieldForceModelContext<T> {

  37.     /** A = sqrt(μ * a). */
  38.     private final T A;

  39.     /** &Chi; = 1 / sqrt(1 - e²) = 1 / B. */
  40.     private final T chi;

  41.     /** &Chi;². */
  42.     private final T chi2;

  43.     // Common factors from equinoctial coefficients
  44.     /** 2 * a / A . */
  45.     private final T ax2oA;

  46.     /** 1 / (A * B) . */
  47.     private final T ooAB;

  48.     /** B / A . */
  49.     private final T BoA;

  50.     /** B / (A * (1 + B)) . */
  51.     private final T BoABpo;

  52.     /** C / (2 * A * B) . */
  53.     private final T Co2AB;

  54.     /** μ / a . */
  55.     private final T muoa;

  56.     /** R / a . */
  57.     private final T roa;

  58.     /** Keplerian mean motion. */
  59.     private final T n;

  60.     /** Direction cosine α. */
  61.     private final T alpha;

  62.     /** Direction cosine β. */
  63.     private final T beta;

  64.     /** Direction cosine γ. */
  65.     private final T gamma;

  66.     /** Transform from body-fixed frame to inertial frame. */
  67.     private final FieldStaticTransform<T> bodyFixedToInertialTransform;

  68.     /**
  69.      * Constructor.
  70.      *
  71.      * @param auxiliaryElements auxiliary elements related to the current orbit
  72.      * @param centralBodyFixedFrame  rotating body frame
  73.      * @param provider   provider for spherical harmonics
  74.      * @param parameters values of the force model parameters
  75.      */
  76.     FieldDSSTGravityContext(final FieldAuxiliaryElements<T> auxiliaryElements,
  77.                             final Frame centralBodyFixedFrame,
  78.                             final UnnormalizedSphericalHarmonicsProvider provider,
  79.                             final T[] parameters) {

  80.         super(auxiliaryElements);

  81.         // µ
  82.         final T mu = parameters[0];

  83.         // Semi-major axis
  84.         final T a = auxiliaryElements.getSma();

  85.         // Keplerian Mean Motion
  86.         final T absA = FastMath.abs(a);
  87.         this.n = FastMath.sqrt(mu.divide(absA)).divide(absA);

  88.         // A = sqrt(µ * |a|)
  89.         this.A = FastMath.sqrt(mu.multiply(absA));

  90.         // &Chi; = 1 / B
  91.         final T B = auxiliaryElements.getB();
  92.         this.chi = auxiliaryElements.getB().reciprocal();
  93.         this.chi2 = chi.multiply(chi);

  94.         // Common factors from equinoctial coefficients
  95.         // 2 * a / A
  96.         this.ax2oA = a.divide(A).multiply(2.);
  97.         // B / A
  98.         this.BoA = B.divide(A);;
  99.         // 1 / AB
  100.         this.ooAB = A.multiply(B).reciprocal();;
  101.         // C / 2AB
  102.         this.Co2AB =  auxiliaryElements.getC().multiply(ooAB).divide(2.);;
  103.         // B / (A * (1 + B))
  104.         this.BoABpo = BoA.divide(B.add(1.));
  105.         // &mu / a
  106.         this.muoa = mu.divide(a);
  107.         // R / a
  108.         this.roa = a.divide(provider.getAe()).reciprocal();


  109.         // If (centralBodyFrame == null), then centralBodyFrame = orbit frame (see DSSTZonal constructors for more on this).
  110.         final Frame internalCentralBodyFrame = centralBodyFixedFrame == null ? auxiliaryElements.getFrame() : centralBodyFixedFrame;

  111.         // Transform from body-fixed frame (typically ITRF) to inertial frame
  112.         bodyFixedToInertialTransform = internalCentralBodyFrame.
  113.                         getStaticTransformTo(auxiliaryElements.getFrame(), auxiliaryElements.getDate());

  114.         final FieldVector3D<T> zB = bodyFixedToInertialTransform.transformVector(Vector3D.PLUS_K);

  115.         // Direction cosines for central body [Eq. 2.1.9-(1)]
  116.         alpha = FieldVector3D.dotProduct(zB, auxiliaryElements.getVectorF());
  117.         beta  = FieldVector3D.dotProduct(zB, auxiliaryElements.getVectorG());
  118.         gamma = FieldVector3D.dotProduct(zB, auxiliaryElements.getVectorW());
  119.     }

  120.     /** A = sqrt(μ * a).
  121.      * @return A
  122.      */
  123.     public T getA() {
  124.         return A;
  125.     }

  126.     /** Get &Chi; = 1 / sqrt(1 - e²) = 1 / B.
  127.      * @return chi
  128.      */
  129.     public T getChi() {
  130.         return chi;
  131.     }

  132.     /** Get &Chi;².
  133.      * @return chi2
  134.      */
  135.     public T getChi2() {
  136.         return chi2;
  137.     }

  138.     /** Getter for the ax2oA.
  139.      * @return the ax2oA
  140.      */
  141.     public T getAx2oA() {
  142.         return ax2oA;
  143.     }

  144.     /** Get ooAB = 1 / (A * B).
  145.      * @return ooAB
  146.      */
  147.     public T getOoAB() {
  148.         return ooAB;
  149.     }

  150.     /** Get B / A.
  151.      * @return BoA
  152.      */
  153.     public T getBoA() {
  154.         return BoA;
  155.     }

  156.     /** Get BoABpo = B / A(1 + B).
  157.      * @return BoABpo
  158.      */
  159.     public T getBoABpo() {
  160.         return BoABpo;
  161.     }

  162.     /** Get Co2AB = C / 2AB.
  163.      * @return Co2AB
  164.      */
  165.     public T getCo2AB() {
  166.         return Co2AB;
  167.     }

  168.     /** Get muoa = μ / a.
  169.      * @return the muoa
  170.      */
  171.     public T getMuoa() {
  172.         return muoa;
  173.     }

  174.     /** Get roa = R / a.
  175.      * @return roa
  176.      */
  177.     public T getRoa() {
  178.         return roa;
  179.     }

  180.     /** Get the Keplerian mean motion.
  181.      * <p>The Keplerian mean motion is computed directly from semi major axis
  182.      * and central acceleration constant.</p>
  183.      * @return Keplerian mean motion in radians per second
  184.      */
  185.     public T getMeanMotion() {
  186.         return n;
  187.     }

  188.     /** Get direction cosine α for central body.
  189.      * @return α
  190.      */
  191.     public T getAlpha() {
  192.         return alpha;
  193.     }

  194.     /** Get direction cosine β for central body.
  195.      * @return β
  196.      */
  197.     public T getBeta() {
  198.         return beta;
  199.     }

  200.     /** Get direction cosine γ for central body.
  201.      * @return the γ
  202.      */
  203.     public T getGamma() {
  204.         return gamma;
  205.     }

  206.     /** Getter for the bodyFixedToInertialTransform.
  207.      * @return the bodyFixedToInertialTransform
  208.      */
  209.     public FieldStaticTransform<T> getBodyFixedToInertialTransform() {
  210.         return bodyFixedToInertialTransform;
  211.     }
  212. }