FieldDSSTThirdBodyDynamicContext.java

  1. /* Copyright 2002-2023 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.util.FastMath;
  21. import org.orekit.bodies.CelestialBody;
  22. import org.orekit.propagation.semianalytical.dsst.utilities.FieldAuxiliaryElements;

  23. /**
  24.  * This class is a container for the common "field" parameters used in {@link DSSTThirdBody}.
  25.  * <p>
  26.  * It performs parameters initialization at each integration step for the third
  27.  * body attraction perturbation. These parameters change for each integration
  28.  * step.
  29.  * </p>
  30.  * @author Bryan Cazabonne
  31.  * @since 12.0
  32.  * @param <T> type of the field elements
  33.  */
  34. public class FieldDSSTThirdBodyDynamicContext<T extends CalculusFieldElement<T>> extends FieldForceModelContext<T> {

  35.     /** Standard gravitational parameter μ for the body in m³/s². */
  36.     private T gm;

  37.     /** Distance from center of mass of the central body to the 3rd body. */
  38.     private T R3;

  39.     /** A = sqrt(μ * a). */
  40.     private T A;

  41.     /** α. */
  42.     private T alpha;

  43.     /** β. */
  44.     private T beta;

  45.     /** γ. */
  46.     private T gamma;

  47.     /** B². */
  48.     private T BB;

  49.     /** B³. */
  50.     private T BBB;

  51.     /** &Chi; = 1 / sqrt(1 - e²) = 1 / B. */
  52.     private T X;

  53.     /** &Chi;². */
  54.     private T XX;

  55.     /** &Chi;³. */
  56.     private T XXX;

  57.     /** -2 * a / A. */
  58.     private T m2aoA;

  59.     /** B / A. */
  60.     private T BoA;

  61.     /** 1 / (A * B). */
  62.     private T ooAB;

  63.     /** -C / (2 * A * B). */
  64.     private T mCo2AB;

  65.     /** B / A(1 + B). */
  66.     private T BoABpo;

  67.     /** mu3 / R3. */
  68.     private T muoR3;

  69.     /** b = 1 / (1 + sqrt(1 - e²)) = 1 / (1 + B).*/
  70.     private T b;

  71.     /** h * &Chi;³. */
  72.     private T hXXX;

  73.     /** k * &Chi;³. */
  74.     private T kXXX;

  75.     /** Keplerian mean motion. */
  76.     private T motion;

  77.     /** Constructor.
  78.      * @param aux auxiliary elements related to the current orbit
  79.      * @param body body the 3rd body to consider
  80.      * @param parameters values of the force model parameters
  81.      */
  82.     public FieldDSSTThirdBodyDynamicContext(final FieldAuxiliaryElements<T> aux,
  83.                                             final CelestialBody body,
  84.                                             final T[] parameters) {

  85.         super(aux);

  86.         // Parameters related to force model drivers
  87.         final T mu = parameters[1];
  88.         A = FastMath.sqrt(mu.multiply(aux.getSma()));
  89.         this.gm = parameters[0];
  90.         final T absA = FastMath.abs(aux.getSma());
  91.         motion = FastMath.sqrt(mu.divide(absA)).divide(absA);

  92.         // Distance from center of mass of the central body to the 3rd body
  93.         final FieldVector3D<T> bodyPos = body.getPVCoordinates(aux.getDate(), aux.getFrame()).getPosition();
  94.         R3 = bodyPos.getNorm();

  95.         // Direction cosines
  96.         final FieldVector3D<T> bodyDir = bodyPos.normalize();
  97.         alpha = bodyDir.dotProduct(aux.getVectorF());
  98.         beta  = bodyDir.dotProduct(aux.getVectorG());
  99.         gamma = bodyDir.dotProduct(aux.getVectorW());

  100.         // &Chi;<sup>-2</sup>.
  101.         BB = aux.getB().multiply(aux.getB());
  102.         // &Chi;<sup>-3</sup>.
  103.         BBB = BB.multiply(aux.getB());

  104.         // b = 1 / (1 + B)
  105.         b = aux.getB().add(1.).reciprocal();

  106.         // &Chi;
  107.         X = aux.getB().reciprocal();
  108.         XX = X.multiply(X);
  109.         XXX = X.multiply(XX);
  110.         // -2 * a / A
  111.         m2aoA = aux.getSma().multiply(-2.).divide(A);
  112.         // B / A
  113.         BoA = aux.getB().divide(A);
  114.         // 1 / AB
  115.         ooAB = (A.multiply(aux.getB())).reciprocal();
  116.         // -C / 2AB
  117.         mCo2AB = aux.getC().multiply(ooAB).divide(2.).negate();
  118.         // B / A(1 + B)
  119.         BoABpo = BoA.divide(aux.getB().add(1.));

  120.         // mu3 / R3
  121.         muoR3 = R3.divide(gm).reciprocal();

  122.         // h * &Chi;³
  123.         hXXX = XXX.multiply(aux.getH());
  124.         // k * &Chi;³
  125.         kXXX = XXX.multiply(aux.getK());

  126.     }

  127.     /** Get A = sqrt(μ * a).
  128.      * @return A
  129.      */
  130.     public T getA() {
  131.         return A;
  132.     }

  133.     /** Get the distance from center of mass of the central body to the 3rd body.
  134.      * @return the distance from center of mass of the central body to the 3rd body
  135.      */
  136.     public T getR3() {
  137.         return R3;
  138.     }

  139.     /** Get direction cosine α for central body.
  140.      * @return α
  141.      */
  142.     public T getAlpha() {
  143.         return alpha;
  144.     }

  145.     /** Get direction cosine β for central body.
  146.      * @return β
  147.      */
  148.     public T getBeta() {
  149.         return beta;
  150.     }

  151.     /** Get direction cosine γ for central body.
  152.      * @return γ
  153.      */
  154.     public T getGamma() {
  155.         return gamma;
  156.     }

  157.     /** Get B².
  158.      * @return B²
  159.      */
  160.     public T getBB() {
  161.         return BB;
  162.     }

  163.     /** Get B³.
  164.      * @return B³
  165.      */
  166.     public T getBBB() {
  167.         return BBB;
  168.     }

  169.     /** Get b = 1 / (1 + sqrt(1 - e²)) = 1 / (1 + B).
  170.      * @return b
  171.      */
  172.     public T getb() {
  173.         return b;
  174.     }

  175.     /** Get &Chi; = 1 / sqrt(1 - e²) = 1 / B.
  176.      * @return &Chi;
  177.      */
  178.     public T getX() {
  179.         return X;
  180.     }

  181.     /** Get &Chi;².
  182.      * @return &Chi;²
  183.      */
  184.     public T getXX() {
  185.         return XX;
  186.     }

  187.     /** Get m2aoA = -2 * a / A.
  188.      * @return m2aoA
  189.      */
  190.     public T getM2aoA() {
  191.         return m2aoA;
  192.     }

  193.     /** Get B / A.
  194.      * @return BoA
  195.      */
  196.     public T getBoA() {
  197.         return BoA;
  198.     }

  199.     /** Get ooAB = 1 / (A * B).
  200.      * @return ooAB
  201.      */
  202.     public T getOoAB() {
  203.         return ooAB;
  204.     }

  205.     /** Get mCo2AB = -C / 2AB.
  206.      * @return mCo2AB
  207.      */
  208.     public T getMCo2AB() {
  209.         return mCo2AB;
  210.     }

  211.     /** Get BoABpo = B / A(1 + B).
  212.      * @return BoABpo
  213.      */
  214.     public T getBoABpo() {
  215.         return BoABpo;
  216.     }

  217.     /** Get muoR3 = mu3 / R3.
  218.      * @return muoR3
  219.      */
  220.     public T getMuoR3() {
  221.         return muoR3;
  222.     }

  223.     /** Get hXXX = h * &Chi;³.
  224.      * @return hXXX
  225.      */
  226.     public T getHXXX() {
  227.         return hXXX;
  228.     }

  229.     /** Get kXXX = h * &Chi;³.
  230.      * @return kXXX
  231.      */
  232.     public T getKXXX() {
  233.         return kXXX;
  234.     }

  235.     /** Get the Keplerian mean motion.
  236.      * <p>The Keplerian mean motion is computed directly from semi major axis
  237.      * and central acceleration constant.</p>
  238.      * @return Keplerian mean motion in radians per second
  239.      */
  240.     public T getMeanMotion() {
  241.         return motion;
  242.     }

  243. }