SP3Segment.java

  1. /* Copyright 2022-2025 Luc Maisonobe
  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.files.sp3;

  18. import org.hipparchus.analysis.interpolation.HermiteInterpolator;
  19. import org.orekit.attitudes.AttitudeProvider;
  20. import org.orekit.attitudes.FrameAlignedProvider;
  21. import org.orekit.files.general.EphemerisFile;
  22. import org.orekit.files.general.EphemerisSegmentPropagator;
  23. import org.orekit.frames.Frame;
  24. import org.orekit.propagation.BoundedPropagator;
  25. import org.orekit.propagation.SpacecraftState;
  26. import org.orekit.time.AbsoluteDate;
  27. import org.orekit.time.ClockModel;
  28. import org.orekit.time.ClockOffset;
  29. import org.orekit.time.SampledClockModel;
  30. import org.orekit.utils.CartesianDerivativesFilter;
  31. import org.orekit.utils.SortedListTrimmer;

  32. import java.util.ArrayList;
  33. import java.util.Collections;
  34. import java.util.List;

  35. /** One segment of an {@link SP3Ephemeris}.
  36.  * @author Thomas Neidhart
  37.  * @author Evan Ward
  38.  * @author Luc Maisonobe
  39.  * @since 12.0
  40.  */
  41. public class SP3Segment implements EphemerisFile.EphemerisSegment<SP3Coordinate> {

  42.     /** Standard gravitational parameter in m³ / s². */
  43.     private final double mu;

  44.     /** Reference frame. */
  45.     private final Frame frame;

  46.     /** Number of points to use for interpolation. */
  47.     private final int interpolationSamples;

  48.     /** Available derivatives. */
  49.     private final CartesianDerivativesFilter filter;

  50.     /** Ephemeris Data. */
  51.     private final List<SP3Coordinate> coordinates;

  52.     /** Simple constructor.
  53.      * @param mu standard gravitational parameter to use for creating
  54.      * {@link org.orekit.orbits.Orbit Orbits} from the ephemeris data.
  55.      * @param frame reference frame
  56.      * @param interpolationSamples number of points to use for interpolation
  57.      * @param filter available derivatives
  58.      */
  59.     public SP3Segment(final double mu, final Frame frame,
  60.                       final int interpolationSamples, final CartesianDerivativesFilter filter) {
  61.         this.mu                   = mu;
  62.         this.frame                = frame;
  63.         this.interpolationSamples = interpolationSamples;
  64.         this.filter               = filter;
  65.         this.coordinates          = new ArrayList<>();
  66.     }

  67.     /** Extract the clock model.
  68.      * <p>
  69.      * If some clock or clock rate are present in the SP3 files as default values (999999.999999), then they
  70.      * filtered out here when building the clock model, so interpolation will work if at least there are
  71.      * some remaining regular values.
  72.      * </p>
  73.      * @return extracted clock model
  74.      * @since 12.1
  75.      */
  76.     public ClockModel extractClockModel() {
  77.         final List<ClockOffset> sample = new ArrayList<>(coordinates.size());
  78.         coordinates.forEach(c -> {
  79.             final AbsoluteDate date   = c.getDate();
  80.             final double       offset = c.getClockCorrection();
  81.             if (!Double.isNaN(offset)) {
  82.                 final double rate = filter.getMaxOrder() > 0 ? c.getClockRateChange() : Double.NaN;
  83.                 sample.add(new ClockOffset(date, offset, rate, Double.NaN));
  84.             }
  85.         });
  86.         return new SampledClockModel(sample, interpolationSamples);
  87.     }

  88.     /** {@inheritDoc} */
  89.     @Override
  90.     public double getMu() {
  91.         return mu;
  92.     }

  93.     /** {@inheritDoc} */
  94.     @Override
  95.     public AbsoluteDate getStart() {
  96.         return coordinates.get(0).getDate();
  97.     }

  98.     /** {@inheritDoc} */
  99.     @Override
  100.     public AbsoluteDate getStop() {
  101.         return coordinates.get(coordinates.size() - 1).getDate();
  102.     }

  103.     /** {@inheritDoc} */
  104.     @Override
  105.     public Frame getFrame() {
  106.         return frame;
  107.     }

  108.     /** {@inheritDoc} */
  109.     @Override
  110.     public int getInterpolationSamples() {
  111.         return interpolationSamples;
  112.     }

  113.     /** {@inheritDoc} */
  114.     @Override
  115.     public CartesianDerivativesFilter getAvailableDerivatives() {
  116.         return filter;
  117.     }

  118.     /** {@inheritDoc} */
  119.     @Override
  120.     public List<SP3Coordinate> getCoordinates() {
  121.         return Collections.unmodifiableList(this.coordinates);
  122.     }

  123.     /** Adds a new P/V coordinate.
  124.      * @param coord the P/V coordinate of the satellite
  125.      */
  126.     public void addCoordinate(final SP3Coordinate coord) {
  127.         coordinates.add(coord);
  128.     }

  129.     /** {@inheritDoc} */
  130.     @Override
  131.     public BoundedPropagator getPropagator() {
  132.         return new PropagatorWithClock(new FrameAlignedProvider(getInertialFrame()));
  133.     }

  134.     /** {@inheritDoc} */
  135.     @Override
  136.     public BoundedPropagator getPropagator(final AttitudeProvider attitudeProvider) {
  137.         return new PropagatorWithClock(attitudeProvider);
  138.     }

  139.     /** Propagator including clock.
  140.      * @since 12.1
  141.      */
  142.     private class PropagatorWithClock extends EphemerisSegmentPropagator<SP3Coordinate> {

  143.         /** Trimmer for coordinates list. */
  144.         private final SortedListTrimmer trimmer;

  145.         /** Simple constructor.
  146.          * @param attitudeProvider attitude porovider
  147.          */
  148.         PropagatorWithClock(final AttitudeProvider attitudeProvider) {
  149.             super(SP3Segment.this, attitudeProvider);
  150.             this.trimmer = new SortedListTrimmer(getInterpolationSamples());
  151.         }

  152.         /** {@inheritDoc} */
  153.         @Override
  154.         public SpacecraftState updateAdditionalData(final SpacecraftState original) {

  155.             final HermiteInterpolator interpolator = new HermiteInterpolator();

  156.             // Fill interpolator with sample
  157.             trimmer.
  158.                 getNeighborsSubList(original.getDate(), coordinates).
  159.                 forEach(c -> {
  160.                     final double deltaT = c.getDate().durationFrom(original.getDate());
  161.                     if (filter.getMaxOrder() < 1) {
  162.                         // we use only clock offset
  163.                         interpolator.addSamplePoint(deltaT,
  164.                                                     new double[] { c.getClockCorrection() });
  165.                     } else {
  166.                         // we use both clock offset and clock rate
  167.                         interpolator.addSamplePoint(deltaT,
  168.                                                     new double[] { c.getClockCorrection() },
  169.                                                     new double[] { c.getClockRateChange() });
  170.                     }
  171.                 });

  172.             // perform interpolation (we get derivatives even if we used only clock offset)
  173.             final double[][] derivatives = interpolator.derivatives(0.0, 1);

  174.             // add the clock offset and its first derivative
  175.             return super.updateAdditionalData(original).
  176.                 addAdditionalData(SP3Utils.CLOCK_ADDITIONAL_STATE, derivatives[0]).
  177.                 addAdditionalStateDerivative(SP3Utils.CLOCK_ADDITIONAL_STATE, derivatives[1]);

  178.         }

  179.     }

  180. }