AberrationModifier.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.estimation.measurements.modifiers;

  18. import java.util.Collections;
  19. import java.util.HashMap;
  20. import java.util.List;
  21. import java.util.Map;

  22. import org.hipparchus.Field;
  23. import org.hipparchus.analysis.differentiation.Gradient;
  24. import org.hipparchus.analysis.differentiation.GradientField;
  25. import org.hipparchus.geometry.euclidean.threed.FieldVector3D;
  26. import org.hipparchus.geometry.euclidean.threed.Vector3D;
  27. import org.hipparchus.util.FastMath;
  28. import org.hipparchus.util.MathUtils;
  29. import org.orekit.annotation.DefaultDataContext;
  30. import org.orekit.data.DataContext;
  31. import org.orekit.errors.OrekitException;
  32. import org.orekit.errors.OrekitMessages;
  33. import org.orekit.estimation.measurements.AngularRaDec;
  34. import org.orekit.estimation.measurements.EstimatedMeasurement;
  35. import org.orekit.estimation.measurements.EstimatedMeasurementBase;
  36. import org.orekit.estimation.measurements.EstimationModifier;
  37. import org.orekit.estimation.measurements.GroundStation;
  38. import org.orekit.frames.FieldTransform;
  39. import org.orekit.frames.Frame;
  40. import org.orekit.time.AbsoluteDate;
  41. import org.orekit.time.FieldAbsoluteDate;
  42. import org.orekit.utils.Constants;
  43. import org.orekit.utils.FieldPVCoordinates;
  44. import org.orekit.utils.PVCoordinates;
  45. import org.orekit.utils.ParameterDriver;
  46. import org.orekit.utils.TimeSpanMap;
  47. import org.orekit.utils.TimeStampedFieldPVCoordinates;


  48. /**
  49.  * Class modifying theoretical angular measurement with (the inverse of) stellar aberration.
  50.  * <p>
  51.  * This class implements equation 3.252-3 from Seidelmann, "Explanatory Supplement to the Astronmical Almanac", 1992.
  52.  *
  53.  * @author Mark Rutten
  54.  */
  55. public class AberrationModifier implements EstimationModifier<AngularRaDec> {

  56.     /** Data context. */
  57.     private final DataContext dataContext;

  58.     /** Empty constructor.
  59.      * <p>
  60.      * This constructor uses the {@link DefaultDataContext default data context}
  61.      * </p>
  62.      * @since 12.0
  63.      * @see #AberrationModifier(DataContext)
  64.      */
  65.     @DefaultDataContext
  66.     public AberrationModifier() {
  67.         this(DataContext.getDefault());
  68.     }

  69.     /** Constructor.
  70.      * @param dataContext data context
  71.      * @since 12.0.1
  72.      */
  73.     public AberrationModifier(final DataContext dataContext) {
  74.         this.dataContext = dataContext;
  75.     }

  76.     /** {@inheritDoc} */
  77.     @Override
  78.     public String getEffectName() {
  79.         return "aberration";
  80.     }

  81.     /** Natural to proper correction for aberration of light.
  82.      * @param naturalRaDec the "natural" direction (in barycentric coordinates)
  83.      * @param station      the observer ground station
  84.      * @param date         the date of the measurement
  85.      * @param frame        the frame of the measurement
  86.      * @return the "proper" direction (station-relative coordinates)
  87.      */
  88.     @DefaultDataContext
  89.     public static double[] naturalToProper(final double[] naturalRaDec, final GroundStation station,
  90.                                            final AbsoluteDate date, final Frame frame) {
  91.         return naturalToProper(naturalRaDec, station, date, frame, DataContext.getDefault());
  92.     }

  93.     /**
  94.      * Natural to proper correction for aberration of light.
  95.      *
  96.      * @param naturalRaDec the "natural" direction (in barycentric coordinates)
  97.      * @param station      the observer ground station
  98.      * @param date         the date of the measurement
  99.      * @param frame        the frame of the measurement
  100.      * @param context      the data context
  101.      * @return the "proper" direction (station-relative coordinates)
  102.      * @since 12.0.1
  103.      */
  104.     public static double[] naturalToProper(final double[] naturalRaDec, final GroundStation station,
  105.                                            final AbsoluteDate date, final Frame frame, final DataContext context) {

  106.         ensureFrameIsPseudoInertial(frame);

  107.         // Velocity of station relative to barycentre (units of c)
  108.         final PVCoordinates baryPV = context.getCelestialBodies().getSolarSystemBarycenter().getPVCoordinates(date, frame);
  109.         final PVCoordinates stationPV = station.getBaseFrame().getPVCoordinates(date, frame);
  110.         final Vector3D stationBaryVel = stationPV.getVelocity()
  111.                 .subtract(baryPV.getVelocity())
  112.                 .scalarMultiply(1.0 / Constants.SPEED_OF_LIGHT);

  113.         // Delegate to private method
  114.         return lorentzVelocitySum(naturalRaDec, stationBaryVel);
  115.     }

  116.     /**
  117.      * Proper to natural correction for aberration of light.
  118.      *
  119.      * @param properRaDec the "proper" direction (station-relative coordinates)
  120.      * @param station     the observer ground station
  121.      * @param date        the date of the measurement
  122.      * @param frame       the frame of the measurement
  123.      * @return the "natural" direction (in barycentric coordinates)
  124.      */
  125.     @DefaultDataContext
  126.     public static double[] properToNatural(final double[] properRaDec, final GroundStation station,
  127.                                            final AbsoluteDate date, final Frame frame) {
  128.         return properToNatural(properRaDec, station, date, frame, DataContext.getDefault());
  129.     }

  130.     /**
  131.      * Proper to natural correction for aberration of light.
  132.      *
  133.      * @param properRaDec the "proper" direction (station-relative coordinates)
  134.      * @param station     the observer ground station
  135.      * @param date        the date of the measurement
  136.      * @param frame       the frame of the measurement
  137.      * @param context     the data context
  138.      * @return the "natural" direction (in barycentric coordinates)
  139.      * @since 12.0.1
  140.      */
  141.     public static double[] properToNatural(final double[] properRaDec, final GroundStation station,
  142.                                            final AbsoluteDate date, final Frame frame, final DataContext context) {

  143.         // Check measurement frame is inertial
  144.         ensureFrameIsPseudoInertial(frame);

  145.         // Velocity of barycentre relative to station (units of c)
  146.         final PVCoordinates baryPV = context.getCelestialBodies().getSolarSystemBarycenter().getPVCoordinates(date, frame);
  147.         final PVCoordinates stationPV = station.getBaseFrame().getPVCoordinates(date, frame);
  148.         final Vector3D baryStationVel = baryPV.getVelocity()
  149.                 .subtract(stationPV.getVelocity())
  150.                 .scalarMultiply(1.0 / Constants.SPEED_OF_LIGHT);

  151.         // Delegate to private method
  152.         return lorentzVelocitySum(properRaDec, baryStationVel);
  153.     }

  154.     /**
  155.      * Relativistic sum of velocities.
  156.      * This is based on equation 3.252-3 from Seidelmann, "Explanatory Supplement to the Astronmical Almanac", 1992.
  157.      *
  158.      * @param raDec    the direction to transform
  159.      * @param velocity the velocity (units of c)
  160.      * @return the transformed direction
  161.      */
  162.     private static double[] lorentzVelocitySum(final double[] raDec, final Vector3D velocity) {

  163.         // Measurement as unit vector
  164.         final Vector3D direction = new Vector3D(raDec[0], raDec[1]);

  165.         // Coefficients for calculations
  166.         final double inverseBeta = FastMath.sqrt(1.0 - velocity.getNormSq());
  167.         final double velocityScale = 1.0 + direction.dotProduct(velocity) / (1.0 + inverseBeta);

  168.         // From Seidelmann, equation 3.252-3 (unnormalised)
  169.         final Vector3D transformDirection = (direction.scalarMultiply(inverseBeta))
  170.                 .add(velocity.scalarMultiply(velocityScale));
  171.         return new double[] {transformDirection.getAlpha(), transformDirection.getDelta()};
  172.     }

  173.     /**
  174.      * Natural to proper correction for aberration of light.
  175.      *
  176.      * @param naturalRaDec      the "natural" direction (in barycentric coordinates)
  177.      * @param stationToInertial the transform from station to inertial coordinates
  178.      * @param frame             the frame of the measurement
  179.      * @return the "proper" direction (station-relative coordinates)
  180.      */
  181.     @DefaultDataContext
  182.     public static Gradient[] fieldNaturalToProper(final Gradient[] naturalRaDec,
  183.                                                   final FieldTransform<Gradient> stationToInertial,
  184.                                                   final Frame frame) {
  185.         return fieldNaturalToProper(naturalRaDec, stationToInertial, frame, DataContext.getDefault());
  186.     }

  187.     /**
  188.      * Natural to proper correction for aberration of light.
  189.      *
  190.      * @param naturalRaDec      the "natural" direction (in barycentric coordinates)
  191.      * @param stationToInertial the transform from station to inertial coordinates
  192.      * @param frame             the frame of the measurement
  193.      * @param context           the data context
  194.      * @return the "proper" direction (station-relative coordinates)
  195.      * @since 12.0.1
  196.      */
  197.     public static Gradient[] fieldNaturalToProper(final Gradient[] naturalRaDec,
  198.                                                   final FieldTransform<Gradient> stationToInertial,
  199.                                                   final Frame frame,
  200.                                                   final DataContext context) {

  201.         // Check measurement frame is inertial
  202.         ensureFrameIsPseudoInertial(frame);

  203.         // Set up field
  204.         final Field<Gradient> field = naturalRaDec[0].getField();
  205.         final FieldVector3D<Gradient> zero = FieldVector3D.getZero(field);
  206.         final FieldAbsoluteDate<Gradient> date = stationToInertial.getFieldDate();

  207.         // Barycentre in inertial coordinates
  208.         final FieldPVCoordinates<Gradient> baryPV = context.getCelestialBodies().getSolarSystemBarycenter().getPVCoordinates(date, frame);

  209.         // Station in inertial coordinates
  210.         final TimeStampedFieldPVCoordinates<Gradient> stationPV =
  211.                 stationToInertial.transformPVCoordinates(new TimeStampedFieldPVCoordinates<>(date, zero, zero, zero));

  212.         // Velocity of station relative to barycentre (units of c)
  213.         final FieldVector3D<Gradient> stationBaryVel = stationPV.getVelocity()
  214.                 .subtract(baryPV.getVelocity())
  215.                 .scalarMultiply(1.0 / Constants.SPEED_OF_LIGHT);

  216.         return fieldLorentzVelocitySum(naturalRaDec, stationBaryVel);
  217.     }

  218.     /**
  219.      * Proper to natural correction for aberration of light.
  220.      *
  221.      * @param properRaDec       the "proper" direction (station-relative coordinates)
  222.      * @param stationToInertial the transform from station to inertial coordinates
  223.      * @param frame             the frame of the measurement
  224.      * @return the "natural" direction (in barycentric coordinates)
  225.      */
  226.     @DefaultDataContext
  227.     public static Gradient[] fieldProperToNatural(final Gradient[] properRaDec,
  228.                                                   final FieldTransform<Gradient> stationToInertial,
  229.                                                   final Frame frame) {
  230.         return fieldProperToNatural(properRaDec, stationToInertial, frame, DataContext.getDefault());
  231.     }

  232.     /**
  233.      * Proper to natural correction for aberration of light.
  234.      *
  235.      * @param properRaDec       the "proper" direction (station-relative coordinates)
  236.      * @param stationToInertial the transform from station to inertial coordinates
  237.      * @param frame             the frame of the measurement
  238.      * @param context           the data context
  239.      * @return the "natural" direction (in barycentric coordinates)
  240.      * @since 12.0.1
  241.      */
  242.     public static Gradient[] fieldProperToNatural(final Gradient[] properRaDec,
  243.                                                   final FieldTransform<Gradient> stationToInertial,
  244.                                                   final Frame frame,
  245.                                                   final DataContext context) {

  246.         // Check measurement frame is inertial
  247.         ensureFrameIsPseudoInertial(frame);

  248.         // Set up field
  249.         final Field<Gradient> field = properRaDec[0].getField();
  250.         final FieldVector3D<Gradient> zero = FieldVector3D.getZero(field);
  251.         final FieldAbsoluteDate<Gradient> date = stationToInertial.getFieldDate();

  252.         // Barycentre in inertial coordinates
  253.         final FieldPVCoordinates<Gradient> baryPV = context.getCelestialBodies().getSolarSystemBarycenter().getPVCoordinates(date, frame);

  254.         // Station in inertial coordinates
  255.         final TimeStampedFieldPVCoordinates<Gradient> stationPV =
  256.                 stationToInertial.transformPVCoordinates(new TimeStampedFieldPVCoordinates<>(date, zero, zero, zero));

  257.         // Velocity of barycentre relative to station (units of c)
  258.         final FieldVector3D<Gradient> stationBaryVel = stationPV.getVelocity()
  259.                 .negate()
  260.                 .add(baryPV.getVelocity())
  261.                 .scalarMultiply(1.0 / Constants.SPEED_OF_LIGHT);

  262.         return fieldLorentzVelocitySum(properRaDec, stationBaryVel);
  263.     }

  264.     /**
  265.      * Relativistic sum of velocities.
  266.      * This is based on equation 3.252-3 from Seidelmann, "Explanatory Supplement to the Astronmical Almanac", 1992.
  267.      *
  268.      * @param raDec    the direction to transform
  269.      * @param velocity the velocity (units of c)
  270.      * @return the transformed direction
  271.      */
  272.     private static Gradient[] fieldLorentzVelocitySum(final Gradient[] raDec,
  273.                                                       final FieldVector3D<Gradient> velocity) {

  274.         // Measurement as unit vector
  275.         final FieldVector3D<Gradient> direction = new FieldVector3D<>(raDec[0], raDec[1]);

  276.         // Coefficients for calculations
  277.         final Gradient inverseBeta = (velocity.getNormSq().negate().add(1.0)).sqrt();
  278.         final Gradient velocityScale = (direction.dotProduct(velocity)).divide(inverseBeta.add(1.0)).add(1.0);

  279.         // From Seidelmann, equation 3.252-3 (unnormalised)
  280.         final FieldVector3D<Gradient> transformDirection = (direction.scalarMultiply(inverseBeta))
  281.                 .add(velocity.scalarMultiply(velocityScale));
  282.         return new Gradient[] {transformDirection.getAlpha(), transformDirection.getDelta()};
  283.     }


  284.     /** {@inheritDoc} */
  285.     @Override
  286.     public List<ParameterDriver> getParametersDrivers() {
  287.         return Collections.emptyList();
  288.     }

  289.     /** {@inheritDoc} */
  290.     @Override
  291.     public void modifyWithoutDerivatives(final EstimatedMeasurementBase<AngularRaDec> estimated) {

  292.         // Observation date
  293.         final AbsoluteDate date = estimated.getDate();

  294.         // Observation station
  295.         final GroundStation station = estimated.getObservedMeasurement().getStation();

  296.         // Observation frame
  297.         final Frame frame = estimated.getObservedMeasurement().getReferenceFrame();

  298.         // Convert measurement to natural direction
  299.         final double[] estimatedRaDec = estimated.getEstimatedValue();
  300.         final double[] naturalRaDec = properToNatural(estimatedRaDec, station, date, frame, dataContext);

  301.         // Normalise RA
  302.         final double[] observed           = estimated.getObservedValue();
  303.         final double   baseRightAscension = naturalRaDec[0];
  304.         final double   twoPiWrap          = MathUtils.normalizeAngle(baseRightAscension, observed[0]) - baseRightAscension;
  305.         final double   rightAscension     = baseRightAscension + twoPiWrap;

  306.         // New estimated values
  307.         estimated.modifyEstimatedValue(this, rightAscension, naturalRaDec[1]);

  308.     }

  309.     /** {@inheritDoc} */
  310.     @Override
  311.     public void modify(final EstimatedMeasurement<AngularRaDec> estimated) {

  312.         // Observation date
  313.         final AbsoluteDate date = estimated.getDate();

  314.         // Observation frame
  315.         final Frame frame = estimated.getObservedMeasurement().getReferenceFrame();

  316.         // Extract RA/Dec parameters (no state derivatives)
  317.         int nbParams = 6;
  318.         final Map<String, Integer> indices = new HashMap<>();
  319.         for (ParameterDriver driver : estimated.getObservedMeasurement().getParametersDrivers()) {
  320.             if (driver.isSelected()) {
  321.                 for (TimeSpanMap.Span<String> span = driver.getNamesSpanMap().getFirstSpan(); span != null; span = span.next()) {
  322.                     if (!indices.containsKey(span.getData())) {
  323.                         indices.put(span.getData(), nbParams++);
  324.                     }
  325.                 }
  326.             }
  327.         }
  328.         final Field<Gradient> field = GradientField.getField(nbParams);

  329.         // Observation location
  330.         final FieldTransform<Gradient> stationToInertial =
  331.                 estimated.getObservedMeasurement().getStation().getOffsetToInertial(frame, date, nbParams, indices);

  332.         // Convert measurement to natural direction
  333.         final double[] estimatedRaDec = estimated.getEstimatedValue();
  334.         final Gradient[] estimatedRaDecDS = new Gradient[] {
  335.                 field.getZero().add(estimatedRaDec[0]),
  336.                 field.getZero().add(estimatedRaDec[1])
  337.         };
  338.         final Gradient[] naturalRaDec = fieldProperToNatural(estimatedRaDecDS, stationToInertial, frame, dataContext);

  339.         // Normalise RA
  340.         final double[] observed = estimated.getObservedValue();
  341.         final Gradient baseRightAscension = naturalRaDec[0];
  342.         final double twoPiWrap = MathUtils.normalizeAngle(baseRightAscension.getReal(),
  343.                 observed[0]) - baseRightAscension.getReal();
  344.         final Gradient rightAscension = baseRightAscension.add(twoPiWrap);

  345.         // New estimated values
  346.         estimated.modifyEstimatedValue(this, rightAscension.getValue(), naturalRaDec[1].getValue());

  347.         // Derivatives (only parameter, no state)
  348.         final double[] raDerivatives = rightAscension.getGradient();
  349.         final double[] decDerivatives = naturalRaDec[1].getGradient();

  350.         for (final ParameterDriver driver : estimated.getObservedMeasurement().getParametersDrivers()) {
  351.             for (TimeSpanMap.Span<String> span = driver.getNamesSpanMap().getFirstSpan(); span != null; span = span.next()) {
  352.                 final Integer index = indices.get(span.getData());
  353.                 if (index != null) {
  354.                     final double[] parameterDerivative = estimated.getParameterDerivatives(driver);
  355.                     parameterDerivative[0] += raDerivatives[index];
  356.                     parameterDerivative[1] += decDerivatives[index];
  357.                     estimated.setParameterDerivatives(driver, span.getStart(), parameterDerivative[0], parameterDerivative[1]);
  358.                 }
  359.             }
  360.         }
  361.     }

  362.     /**
  363.      * Check that given frame is pseudo-inertial. Throws an error otherwise.
  364.      *
  365.      * @param frame to check
  366.      *
  367.      * @throws OrekitException if given frame is not pseudo-inertial
  368.      */
  369.     private static void ensureFrameIsPseudoInertial(final Frame frame) {
  370.         // Check measurement frame is inertial
  371.         if (!frame.isPseudoInertial()) {
  372.             throw new OrekitException(OrekitMessages.NON_PSEUDO_INERTIAL_FRAME, frame.getName());
  373.         }
  374.     }

  375. }