LeastSquaresTleGenerationAlgorithm.java

  1. /* Copyright 2002-2025 Mark Rutten
  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.  * Mark Rutten 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.analytical.tle.generation;

  18. import org.hipparchus.CalculusFieldElement;
  19. import org.hipparchus.optim.nonlinear.vector.leastsquares.LeastSquaresOptimizer;
  20. import org.hipparchus.optim.nonlinear.vector.leastsquares.LevenbergMarquardtOptimizer;
  21. import org.orekit.annotation.DefaultDataContext;
  22. import org.orekit.data.DataContext;
  23. import org.orekit.frames.Frame;
  24. import org.orekit.orbits.KeplerianOrbit;
  25. import org.orekit.orbits.OrbitType;
  26. import org.orekit.propagation.FieldSpacecraftState;
  27. import org.orekit.propagation.SpacecraftState;
  28. import org.orekit.propagation.analytical.tle.FieldTLE;
  29. import org.orekit.propagation.analytical.tle.TLE;
  30. import org.orekit.propagation.conversion.osc2mean.LeastSquaresConverter;
  31. import org.orekit.propagation.conversion.osc2mean.TLETheory;
  32. import org.orekit.time.TimeScale;
  33. import org.orekit.utils.ParameterDriver;

  34. /**
  35.  * Least squares method to generate a usable TLE from a spacecraft state.
  36.  *
  37.  * @author Mark Rutten
  38.  * @since 12.0
  39.  */
  40. public class LeastSquaresTleGenerationAlgorithm implements TleGenerationAlgorithm {

  41.     /** Default value for maximum number of iterations.*/
  42.     public static final int DEFAULT_MAX_ITERATIONS = 1000;

  43.     /** Osculating to mean orbit converter. */
  44.     private final LeastSquaresConverter converter;

  45.     /** UTC time scale. */
  46.     private final TimeScale utc;

  47.     /**
  48.      * Default constructor.
  49.      * <p>Uses:
  50.      * <ul>
  51.      * <li>the {@link DataContext#getDefault() default data context}</li>
  52.      * <li>{@link #DEFAULT_MAX_ITERATIONS}</li>
  53.      * <li>the {@link LevenbergMarquardtOptimizer}</li>
  54.      * </ul>
  55.      */
  56.     @DefaultDataContext
  57.     public LeastSquaresTleGenerationAlgorithm() {
  58.         this(DEFAULT_MAX_ITERATIONS);
  59.     }

  60.     /**
  61.      * Default constructor.
  62.      * <p>Uses:
  63.      * <ul>
  64.      * <li>the {@link DataContext#getDefault() default data context}</li>
  65.      * <li>the {@link LevenbergMarquardtOptimizer}</li>
  66.      * </ul>
  67.      * @param maxIterations maximum number of iterations for convergence
  68.      */
  69.     @DefaultDataContext
  70.     public LeastSquaresTleGenerationAlgorithm(final int maxIterations) {
  71.         this(maxIterations, DataContext.getDefault().getTimeScales().getUTC(),
  72.              DataContext.getDefault().getFrames().getTEME());
  73.     }

  74.     /**
  75.      * Constructor.
  76.      * <p>Uses the {@link LevenbergMarquardtOptimizer}.</p>
  77.      * @param maxIterations maximum number of iterations for convergence
  78.      * @param utc  UTC time scale
  79.      * @param teme TEME frame
  80.      */
  81.     public LeastSquaresTleGenerationAlgorithm(final int maxIterations,
  82.                                               final TimeScale utc,
  83.                                               final Frame teme) {
  84.         this(utc, teme, new LeastSquaresConverter(new TLETheory(utc, teme),
  85.                                                   new LevenbergMarquardtOptimizer(),
  86.                                                   LeastSquaresConverter.DEFAULT_THRESHOLD,
  87.                                                   maxIterations));
  88.     }

  89.     /**
  90.      * Constructor.
  91.      * <p>Enables to select the {@link LeastSquaresOptimizer optimizer}
  92.      * for the {@link LeastSquaresConverter least-squares converter}.</p>
  93.      * @param utc  UTC time scale
  94.      * @param teme TEME frame
  95.      * @param converter osculating to mean orbit converter using a least-squares algorithm
  96.      */
  97.     public LeastSquaresTleGenerationAlgorithm(final TimeScale utc,
  98.                                               final Frame teme,
  99.                                               final LeastSquaresConverter converter) {
  100.         this.utc       = utc;
  101.         this.converter = converter;
  102.         converter.setMeanTheory(new TLETheory(utc, teme));
  103.     }

  104.     /** {@inheritDoc} */
  105.     @Override
  106.     public TLE generate(final SpacecraftState state, final TLE templateTLE) {
  107.         final KeplerianOrbit mean = (KeplerianOrbit) OrbitType.KEPLERIAN.convertType(converter.convertToMean(state.getOrbit()));
  108.         final TLE tle = TleGenerationUtil.newTLE(mean, templateTLE, templateTLE.getBStar(mean.getDate()), utc);
  109.         // reset estimated parameters from template to generated tle
  110.         for (final ParameterDriver templateDrivers : templateTLE.getParametersDrivers()) {
  111.             if (templateDrivers.isSelected()) {
  112.                 // set to selected for the new TLE
  113.                 tle.getParameterDriver(templateDrivers.getName()).setSelected(true);
  114.             }
  115.         }
  116.         return tle;
  117.     }

  118.     /**
  119.      * Get the Root Mean Square of the TLE estimation.
  120.      * <p>
  121.      * Be careful that the RMS is updated each time the
  122.      * {@link LeastSquaresTleGenerationAlgorithm#generate(SpacecraftState, TLE)}
  123.      * method is called.
  124.      * </p>
  125.      * @return the RMS
  126.      */
  127.     public double getRms() {
  128.         return converter.getRMS();
  129.     }

  130.     /** {@inheritDoc} */
  131.     @Override
  132.     public <T extends CalculusFieldElement<T>> FieldTLE<T> generate(final FieldSpacecraftState<T> state,
  133.                                                                     final FieldTLE<T> templateTLE) {
  134.         throw new UnsupportedOperationException();
  135.     }
  136. }