AberrationModifier.java

  1. /* Copyright 2002-2024 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.     /**
  70.      * Constructor.
  71.      * @param dataContext data context
  72.      * @since 12.0.1
  73.      */
  74.     public AberrationModifier(final DataContext dataContext) {
  75.         this.dataContext = dataContext;
  76.     }

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

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

  104.         ensureFrameIsPseudoInertial(frame);

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

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

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

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

  141.         // Check measurement frame is inertial
  142.         ensureFrameIsPseudoInertial(frame);

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

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

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

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

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

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

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

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

  199.         // Check measurement frame is inertial
  200.         ensureFrameIsPseudoInertial(frame);

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

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

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

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

  214.         return fieldLorentzVelocitySum(naturalRaDec, stationBaryVel);
  215.     }

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

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

  244.         // Check measurement frame is inertial
  245.         ensureFrameIsPseudoInertial(frame);

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

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

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

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

  260.         return fieldLorentzVelocitySum(properRaDec, stationBaryVel);
  261.     }

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

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

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

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


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

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

  290.         // Observation date
  291.         final AbsoluteDate date = estimated.getDate();

  292.         // Observation station
  293.         final GroundStation station = estimated.getObservedMeasurement().getStation();

  294.         // Observation frame
  295.         final Frame frame = estimated.getObservedMeasurement().getReferenceFrame();

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

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

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

  306.     }

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

  310.         // Observation date
  311.         final AbsoluteDate date = estimated.getDate();

  312.         // Observation frame
  313.         final Frame frame = estimated.getObservedMeasurement().getReferenceFrame();

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

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

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

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

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

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

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

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

  373. }