observe.h#
Header file for the observable computations.
- Author
- Rahil Makadia makadia2@illinois.edu 
LICENSE#
Copyright (C) 2022-2025 Rahil Makadia
This program is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 3 of the License, or (at your option) any later version.
This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
You should have received a copy of the GNU General Public License along with this program; if not, see https://www.gnu.org/licenses.
Functions
- 
void get_glb_correction(PropSimulation *propSim, const size_t &interpIdx, const real &tInterpGeom, std::vector<real> &xInterpApparentBary)#
- Compute the correction to the apparent state of the body due to the gravitational light bending. - Parameters:
- propSim – [in] PropSimulation object for the integration. 
- interpIdx – [in] Index of the next interpolation time. 
- tInterpGeom – [in] Time to interpolate to. 
- xInterpApparentBary – [out] Apparent state vector of the target body. 
 
 
- 
void get_measurement(PropSimulation *propSim, const size_t &interpIdx, const real tInterpGeom, const std::vector<real> &xInterpGeom, const std::vector<real> &xInterpApparent)#
- Get the relevant measurement (optical/radar) for a given measurement time. - Parameters:
- propSim – [in] PropSimulation object for the integration. 
- interpIdx – [in] Index of the next interpolation time. 
- tInterpGeom – [in] Time to interpolate to. 
- xInterpGeom – [in] Geometric state vector of the target body at the interpolation time. 
- xInterpApparent – [in] Apparent state vector of the target body at the interpolation time. 
 
 
- 
void get_optical_measurement(PropSimulation *propSim, const std::vector<real> &xInterpApparent, std::vector<real> &opticalMeasurement, std::vector<real> &opticalMeasurementDot, std::vector<real> &opticalPartials)#
- Get the optical measurement and partials. - Parameters:
- propSim – [in] PropSimulation object for the integration. 
- xInterpApparent – [in] Apparent state vector of the target body. 
- opticalMeasurement – [out] Optical measurement (RA and Dec). 
- opticalMeasurementDot – [out] Time derivative of the optical measurement. 
- opticalPartials – [out] Partials of the optical measurement. 
 
 
- 
void get_photocenter_correction(PropSimulation *propSim, const size_t &interpIdx, const real &tInterpGeom, const std::vector<real> &xInterpApparent, std::vector<real> &photocenterCorr)#
- Get the photocenter-barycenter correction for an optical measurement. - Parameters:
- propSim – [in] PropSimulation object for the integration. 
- interpIdx – [in] Index of the next interpolation time. 
- tInterpGeom – [in] Time to interpolate to. 
- xInterpApparent – [in] Apparent state vector of the target body. 
- photocenterCorr – [out] Photocenter-barycenter correction to the optical measurement. 
 
 
- 
void get_radar_measurement(PropSimulation *propSim, const size_t &interpIdx, const real tInterpGeom, const std::vector<real> &xInterpGeom, std::vector<real> &radarMeasurement, std::vector<real> &radarPartials)#
- Get the radar measurement and partials. - Parameters:
- propSim – [in] PropSimulation object for the integration. 
- interpIdx – [in] Index of the next interpolation time. 
- tInterpGeom – [in] Time to interpolate to. 
- xInterpGeom – [in] Geometric state vector of the target body at the interpolation time. 
- radarMeasurement – [out] Radar measurement (range/Doppler). 
- radarPartials – [out] Partials of the radar measurement. 
 
 
- 
void get_delay_measurement(PropSimulation *propSim, const size_t &interpIdx, const size_t &i, const real tInterpGeom, const std::vector<real> &xInterpGeom, const real &receiveTimeTDB, real &transmitTimeTDB, std::vector<real> &xObsBaryRcv, std::vector<real> &xTrgtBaryBounce, std::vector<real> &xObsBaryTx, real &delayMeasurement, std::vector<real> &delayPartials)#
- Get the radar delay measurement and partials. - Parameters:
- propSim – [in] PropSimulation object for the integration. 
- interpIdx – [in] Index of the next interpolation time. 
- i – [in] Index of the target body. 
- tInterpGeom – [in] Time to interpolate to. 
- xInterpGeom – [in] Geometric state vector of the target body at the interpolation time. 
- receiveTimeTDB – [in] Time of reception of the radar signal (TDB). 
- transmitTimeTDB – [out] Time of transmission of the radar signal (TDB). 
- xObsBaryRcv – [out] Barycentric state vector of the observer at the reception time. 
- xTrgtBaryBounce – [out] Barycentric state vector of the target body at the bounce time. 
- xObsBaryTx – [out] Barycentric state vector of the observer at the transmission time. 
- delayMeasurement – [out] Radar delay measurement. 
- delayPartials – [out] Partials of the radar delay measurement. 
 
 
- 
void get_delta_delay_relativistic(PropSimulation *propSim, const real &tForSpice, const std::vector<real> &targetState, real &deltaDelayRelativistic)#
- Get the relativistic delay measurement correction. - Parameters:
- propSim – [in] PropSimulation object for the integration. 
- tForSpice – [in] Time for querying SPICE Ephemeris for the Sun and Earth positions. 
- targetState – [in] State vector of the target body. 
- deltaDelayRelativistic – [out] Relativistic delay correction (Shapiro delay). 
 
 
- 
void get_doppler_measurement(PropSimulation *propSim, const size_t &i, const real receiveTimeTDB, const real transmitTimeTDB, const std::vector<real> xObsBaryRcv, const std::vector<real> xTrgtBaryBounce, const std::vector<real> xObsBaryTx, const real transmitFreq, real &dopplerMeasurement, std::vector<real> &dopplerPartials)#
- Get the Doppler measurement and partials. - Parameters:
- propSim – [in] PropSimulation object for the integration. 
- i – [in] Index of the target body. 
- receiveTimeTDB – [in] Time of reception of the radar signal (TDB). 
- transmitTimeTDB – [in] Time of transmission of the radar signal (TDB). 
- xObsBaryRcv – [in] Barycentric state vector of the observer at the reception time. 
- xTrgtBaryBounce – [in] Barycentric state vector of the target body at the bounce time. 
- xObsBaryTx – [in] Barycentric state vector of the observer at the transmission time. 
- transmitFreq – [in] Frequency of the transmitted radar signal. 
- dopplerMeasurement – [out] Doppler measurement. 
- dopplerPartials – [out] Partials of the Doppler measurement. 
 
 
- 
void evaluate_one_interpolation(const PropSimulation *propSim, const real &tInterp, std::vector<real> &xInterp)#
- Interpolate the integrator state for one evaluation time. - Parameters:
- propSim – [in] PropSimulation object for the integration. 
- tInterp – [in] Time to interpolate to. 
- xInterp – [out] Interpolated state vector.