USGSCSM Astro SAR Sensor Model Class

class UsgsAstroSarSensorModel : public csm::RasterGM, public virtual csm::SettableEllipsoid
#include <UsgsAstroSarSensorModel.h>

Public Types

enum LookDirection

Values:

enumerator LEFT
enumerator RIGHT

Public Functions

UsgsAstroSarSensorModel()

Constructor for UsgsAstroSarSensorModel.

@description Initializes a new instance of the UsgsAstroSarSensorModel class, setting up the initial state and logging the creation of the model.

inline ~UsgsAstroSarSensorModel()
void reset()

Resets the sensor model to its default state.

@description Resets all member variables of the UsgsAstroSarSensorModel to their default values. This includes setting identifiers to “Unknown”, setting numerical values to 0 or appropriate defaults, clearing vectors, and preparing the model for a fresh state configuration.

virtual void replaceModelState(const std::string &argState)

Replaces the current model state with a new one.

@description This function takes a JSON string representation of a sensor model’s state, parses it, and updates the current model state accordingly. It resets the model before applying the new state to ensure a clean update. Throws an error if the look direction is unrecognized.

Parameters:

argState – The new model state in JSON string format.

Throws:

csm::Error – If the look direction from the state is invalid.

virtual std::string getModelState() const

Retrieves the current model state as a JSON string.

@description Serializes the current state of the sensor model into a JSON string. This includes all relevant parameters and configurations that define the sensor model’s state. If the look direction is not correctly defined in the model, it throws an error.

Throws:

csm::Error – If the look direction is not properly set in the model.

Returns:

A JSON string representing the current state of the model.

std::string constructStateFromIsd(const std::string imageSupportData, csm::WarningList *list)

Constructs the sensor model state from ISD (Image Support Data) for SAR sensors.

@description Parses the ISD to extract necessary information for constructing the state of a SAR sensor model. This includes setting up the sensor, platform, and image identifiers, geometric properties, and other SAR-specific parameters. The function also handles linear interpolation of instrument and sun states (positions and velocities) based on provided ephemeris times. Warnings are generated for missing or invalid data.

Parameters:
  • imageSupportData – The ISD in string format containing sensor and image acquisition details.

  • warnings – Pointer to a list for recording any warnings encountered during the state construction.

Returns:

A string representing the constructed state in JSON format.

virtual csm::ImageCoord groundToImage(const csm::EcefCoord &groundPt, double desiredPrecision = 0.001, double *achievedPrecision = NULL, csm::WarningList *warnings = NULL) const

Computes the image coordinates (line and sample) for a given ground point using the current model state without adjustments.

@description Calculates the image coordinates that correspond to a given ground point.

Parameters:
  • groundPt – The ground point in ECEF coordinates.

  • desiredPrecision – The desired precision for the ground-to-image calculation.

  • achievedPrecision – Pointer to a double to store the achieved precision.

  • warnings – Pointer to a warning list to capture any warnings.

Returns:

The computed image coordinates.

virtual csm::ImageCoord groundToImage(const csm::EcefCoord &groundPt, const std::vector<double> adjustments, double desired_precision = 0.001, double *achieved_precision = NULL, csm::WarningList *warnings = NULL) const

Computes the image coordinates (line and sample) for a given ground point using the current model state with adjustments.

@description Calculates the image coordinates that correspond to a given ground point, taking into account any adjustments specified. This function finds the time of the closest approach to the ground point and the corresponding slant range, then calculates the ground range using a polynomial defined by the range coefficient set. It accounts for doppler shift to accurately determine the position in the image.

Parameters:
  • groundPt – The ground point in ECEF coordinates.

  • adj – Adjustments to be applied to the sensor model state.

  • desiredPrecision – The desired precision for the ground-to-image calculation.

  • achievedPrecision – Pointer to a double to store the achieved precision.

  • warnings – Pointer to a warning list to capture any warnings.

Throws:

csm::Error – If the computation fails for any reason.

Returns:

The computed image coordinates.

virtual csm::ImageCoordCovar groundToImage(const csm::EcefCoordCovar &groundPt, double desiredPrecision = 0.001, double *achievedPrecision = NULL, csm::WarningList *warnings = NULL) const

Converts a ground point with covariance to an image coordinate with covariance.

@description This function converts a ground point, specified in ECEF coordinates with associated covariance, into an image coordinate with covariance. It utilizes the groundToImage function for the conversion and then propagates the uncertainty from the ground point through the transformation to compute the covariance in the image space.

Parameters:
  • groundPt – A ground point with covariance in ECEF coordinates.

  • desiredPrecision – The desired precision for the conversion.

  • achievedPrecision – Pointer to store the achieved precision.

  • warnings – Warning list to capture any issues during the process.

Returns:

An image coordinate with covariance resulting from the conversion.

virtual csm::EcefCoord imageToGround(const csm::ImageCoord &imagePt, double height, double desiredPrecision = 0.001, double *achievedPrecision = NULL, csm::WarningList *warnings = NULL) const

Converts an image coordinate to a ground point in ECEF coordinates.

@description This function calculates the ground point corresponding to a given image coordinate and height above the ellipsoid. The process involves determining the time of the image acquisition and the ground range, then converting the ground range to slant range, and finally solving for the ground point position using the spacecraft position and velocity.

Parameters:
  • imagePt – The image coordinate (line, sample) for which to find the corresponding ground point.

  • height – The height above the ellipsoid at which to find the ground point.

  • desiredPrecision – The desired precision for the conversion.

  • achievedPrecision – Pointer to store the achieved precision of the conversion.

  • warnings – Warning list to capture any issues during the process.

Returns:

The corresponding ground point in ECEF coordinates.

virtual csm::EcefCoordCovar imageToGround(const csm::ImageCoordCovar &imagePt, double height, double heightVariance, double desiredPrecision = 0.001, double *achievedPrecision = NULL, csm::WarningList *warnings = NULL) const

Converts an image coordinate with covariance to a ground point with covariance.

@description Converts an image coordinate with covariance to a ground point in ECEF coordinates, considering the height and its variance, sensor modeling, and unmodeled errors for covariance computation.

Parameters:
  • imagePt – Image coordinate with covariance.

  • height – Height above the ellipsoid for the ground point.

  • heightVariance – Variance of the height above the ellipsoid.

  • desiredPrecision – Desired precision for the conversion.

  • achievedPrecision – Pointer to store the achieved precision.

  • warnings – Warning list for issues during the process.

Returns:

Ground point in ECEF coordinates with covariance.

virtual csm::EcefLocus imageToProximateImagingLocus(const csm::ImageCoord &imagePt, const csm::EcefCoord &groundPt, double desiredPrecision = 0.001, double *achievedPrecision = NULL, csm::WarningList *warnings = NULL) const

Computes proximate imaging locus for an image point and ground point.

@description Calculates the closest point on the imaging locus to a ground point, defining the imaging locus as the intersection of the Doppler cone with the Earth’s surface.

Parameters:
  • imagePt – Image coordinate (line, sample).

  • groundPt – Ground point in ECEF coordinates.

  • desiredPrecision – Desired precision for the conversion.

  • achievedPrecision – Pointer to store the achieved precision.

  • warnings – Warning list for issues during the process.

Returns:

Locus of points forming the proximate imaging locus, including the direction vector.

virtual csm::EcefLocus imageToRemoteImagingLocus(const csm::ImageCoord &imagePt, double desiredPrecision = 0.001, double *achievedPrecision = NULL, csm::WarningList *warnings = NULL) const

Computes the remote imaging locus for an image point.

@description Calculates the imaging locus for an image coordinate by determining the direction from the sensor position through the image point towards the ground.

Parameters:
  • imagePt – Image coordinate (line, sample) for locus calculation.

  • desiredPrecision – Desired precision for the calculation.

  • achievedPrecision – Pointer to store achieved precision.

  • warnings – Warning list for issues during the process.

Returns:

Locus of points forming the remote imaging locus, including the direction vector.

virtual csm::ImageCoord getImageStart() const

Returns the starting image coordinate.

@description Always returns the origin point (0.0, 0.0) as the start of the image coordinates.

Returns:

The starting image coordinate.

virtual csm::ImageVector getImageSize() const

Returns the size of the image.

@description Provides the dimensions of the image in terms of the number of lines and samples.

Returns:

Image vector indicating the total number of lines and samples.

virtual std::pair<csm::ImageCoord, csm::ImageCoord> getValidImageRange() const

Provides the valid range of image coordinates.

@description Computes the range of valid image coordinates starting from the origin to the end of the image based on its size.

Returns:

Pair of image coordinates indicating the start and end of the valid image range.

virtual std::pair<double, double> getValidHeightRange() const

Gets the valid range of heights for the model.

@description Returns the minimum and maximum valid elevation values for the sensor model.

Returns:

Pair of doubles indicating the minimum and maximum valid heights.

virtual csm::EcefVector getIlluminationDirection(const csm::EcefCoord &groundPt) const

Computes the illumination direction for a ground point.

@description Determines the unit vector from a ground point towards the sun, indicating the direction of illumination.

Parameters:

groundPt – Ground point in ECEF coordinates.

Returns:

Normalized vector pointing towards the sun from the ground point.

virtual double getImageTime(const csm::ImageCoord &imagePt) const

Computes the image time for a given image coordinate.

@description Calculates the time at which a given line in the image was acquired, based on the starting ephemeris time and exposure duration.

Parameters:

imagePt – Image coordinate (line, sample).

Returns:

The time at which the image line was acquired.

virtual csm::EcefVector getSpacecraftPosition(double time) const

Returns the spacecraft position at a given time without adjustments.

@description Retrieves the ECEF position of the spacecraft at a specific time, without considering any model adjustments.

Parameters:

time – The time at which to retrieve the spacecraft position.

Returns:

Vector representing the spacecraft’s position in ECEF coordinates.

virtual csm::EcefVector getAdjustedSpacecraftPosition(double time, std::vector<double> adj) const

Returns the adjusted spacecraft position at a given time.

@description Retrieves the ECEF position of the spacecraft at a specific time, considering the provided model adjustments.

Parameters:
  • time – The time at which to retrieve the spacecraft position.

  • adj – Vector of adjustments to apply to the model.

Returns:

Vector representing the adjusted spacecraft’s position in ECEF coordinates.

virtual csm::EcefCoord getSensorPosition(const csm::ImageCoord &imagePt) const

Retrieves the sensor position for a given image coordinate.

@description Calculates the sensor position in ECEF coordinates corresponding to a specific image coordinate by first determining the time of the image point and then finding the sensor position at that time.

Parameters:

imagePt – The image coordinate for which to find the sensor position.

Returns:

The ECEF coordinates of the sensor position for the given image coordinate.

virtual csm::EcefCoord getSensorPosition(double time) const

Retrieves the sensor position at a given time.

@description Retrieves the ECEF position of the sensor at a specific time, without applying any model adjustments.

Parameters:

time – The time at which to retrieve the sensor position.

Returns:

The ECEF coordinates of the sensor position at the given time.

virtual csm::EcefCoord getAdjustedSensorPosition(double time, std::vector<double> adjustments) const

Retrieves the adjusted sensor position at a given time.

@description Retrieves the ECEF position of the sensor at a specific time, considering the provided model adjustments.

Parameters:
  • time – The time at which to retrieve the sensor position.

  • adj – Vector of adjustments to apply to the model.

Returns:

The adjusted ECEF coordinates of the sensor position at the given time.

virtual csm::EcefVector getSensorVelocity(const csm::ImageCoord &imagePt) const

Retrieves the sensor velocity for a given image coordinate.

@description Calculates the sensor velocity in ECEF coordinates corresponding to a specific image coordinate by first determining the time of the image point and then finding the sensor velocity at that time. No model adjustments are applied.

Parameters:

imagePt – The image coordinate for which to find the sensor velocity.

Returns:

The ECEF vector representing the sensor velocity for the given image coordinate.

virtual csm::EcefVector getSensorVelocity(double time) const

Retrieves the sensor velocity at a given time.

@description Retrieves the ECEF velocity of the sensor at a specific time, without applying any model adjustments.

Parameters:

time – The time at which to retrieve the sensor velocity.

Returns:

The ECEF vector of the sensor velocity at the given time.

virtual csm::EcefVector getAdjustedSensorVelocity(double time, std::vector<double> adjustments) const

Retrieves the adjusted sensor velocity at a given time.

@description Retrieves the ECEF velocity of the sensor at a specific time, considering the provided model adjustments.

Parameters:
  • time – The time at which to retrieve the sensor velocity.

  • adj – Vector of adjustments to apply to the model.

Returns:

The adjusted ECEF vector of the sensor velocity at the given time.

virtual csm::RasterGM::SensorPartials computeSensorPartials(int index, const csm::EcefCoord &groundPt, double desiredPrecision = 0.001, double *achievedPrecision = NULL, csm::WarningList *warnings = NULL) const

Computes the sensor model partial derivatives for a ground point.

@description Computes the partial derivatives of the image coordinates with respect to the ground point coordinates, using a given ground point in ECEF coordinates, desired precision for the computation, and potential adjustments applied to the sensor model.

Parameters:
  • index – Index of the model parameter for which the partial is being computed.

  • groundPt – Ground point in ECEF coordinates.

  • desiredPrecision – The precision desired in the computation of the partials.

  • achievedPrecision – The precision achieved in the computation.

  • warnings – List to accumulate any warnings.

Returns:

SensorPartials object containing the partial derivatives of the image coordinates with respect to the model parameter.

virtual csm::RasterGM::SensorPartials computeSensorPartials(int index, const csm::ImageCoord &imagePt, const csm::EcefCoord &groundPt, double desiredPrecision = 0.001, double *achievedPrecision = NULL, csm::WarningList *warnings = NULL) const

Computes the sensor model partial derivatives for an image point.

@description Computes the partial derivatives of the image coordinates with respect to a specific sensor model parameter, given an image point, the corresponding ground point, the desired precision for the computation, and any adjustments applied to the sensor model.

Parameters:
  • index – The index of the model parameter for which the partial is being computed.

  • imagePt – The image coordinates (line and sample).

  • groundPt – The ground point in ECEF coordinates associated with the image point.

  • desiredPrecision – The precision desired in the computation of the partials.

  • achievedPrecision – Pointer to a variable where the achieved precision of the computation will be stored. May be NULL if not needed.

  • warnings – List to accumulate any warnings that occur during the computation.

Returns:

SensorPartials object containing the partial derivatives of the image coordinates (line and sample) with respect to the specified model parameter.

virtual std::vector<double> computeGroundPartials(const csm::EcefCoord &groundPt) const

Computes the ground partial derivatives for a ground point.

@description Computes the partial derivatives of the image coordinates with respect to the ground point coordinates, using only the ground point coordinates.

Parameters:

groundPt – Ground point in ECEF coordinates.

Returns:

Vector containing the six partial derivatives of the image coordinates (line and sample) with respect to the ground coordinates (X, Y, Z).

virtual const csm::CorrelationModel &getCorrelationModel() const

Returns the correlation model used by the sensor model.

@description Provides access to the correlation model used for error propagation. SAR sensor models do not model correlation between image points, so this returns a no-correlation model.

Returns:

A reference to a CorrelationModel object representing the correlation model used by the sensor model.

virtual std::vector<double> getUnmodeledCrossCovariance(const csm::ImageCoord &pt1, const csm::ImageCoord &pt2) const

Computes the unmodeled cross covariance between two image points.

@description Computes the covariance matrix between two points in an image due to unmodeled errors. Since unmodeled errors are not considered in SAR sensor models, this function returns a zero matrix.

Parameters:
  • pt1 – First image coordinate.

  • pt2 – Second image coordinate.

Returns:

Vector representing a 2x2 covariance matrix (flattened) filled with zeros.

virtual csm::EcefCoord getReferencePoint() const

Retrieves the current reference point of the model.

@description Returns the reference point for the sensor model in ECEF coordinates. The reference point is typically the ground point closest to the sensor model’s center of projection.

Returns:

The ECEF coordinates of the reference point.

virtual void setReferencePoint(const csm::EcefCoord &groundPt)

Sets the reference point for the sensor model.

@description Updates the sensor model’s reference point using the provided ECEF coordinates.

Parameters:

groundPt – The new reference point in ECEF coordinates.

virtual int getNumParameters() const

Retrieves the number of parameters used by the sensor model.

@description Returns the total number of parameters in the sensor model.

Returns:

The number of parameters.

virtual std::string getParameterName(int index) const

Retrieves the name of a parameter by its index.

@description Provides the name of the parameter corresponding to the given index.

Parameters:

index – The index of the parameter.

Returns:

The name of the parameter at the given index.

virtual std::string getParameterUnits(int index) const

Retrieves the units of a parameter by its index.

@description Provides the units of the parameter corresponding to the given index. For SAR models, all parameters are in meters.

Parameters:

index – The index of the parameter.

Returns:

The units of the parameter at the given index.

virtual bool hasShareableParameters() const

Indicates whether the model parameters are shareable.

@description Determines if the model parameters can be shared between different sensor model instances. SAR sensor models do not support parameter sharing.

Returns:

False, indicating parameters are not shareable.

virtual bool isParameterShareable(int index) const

Indicates if a specific parameter is shareable.

@description Determines if the specific model parameter can be shared between different sensor model instances. SAR sensor models do not support parameter sharing.

Parameters:

index – The index of the parameter.

Returns:

False, indicating the parameter is not shareable.

virtual csm::SharingCriteria getParameterSharingCriteria(int index) const

Retrieves the sharing criteria for a parameter.

@description Provides the criteria under which a parameter can be shared between different sensor model instances. SAR sensor models do not support parameter sharing, so this returns a default SharingCriteria object.

Parameters:

index – The index of the parameter.

Returns:

A SharingCriteria object.

virtual double getParameterValue(int index) const

Retrieves the value of a parameter by its index.

@description Provides the value of the parameter corresponding to the given index.

Parameters:

index – The index of the parameter.

Returns:

The value of the parameter at the given index.

virtual void setParameterValue(int index, double value)

Sets the value of a parameter by its index.

@description Updates the value of the parameter corresponding to the given index.

Parameters:
  • index – The index of the parameter.

  • value – The new value for the parameter.

virtual csm::param::Type getParameterType(int index) const

Retrieves the type of a parameter by its index.

@description Provides the type of the parameter corresponding to the given index.

Parameters:

index – The index of the parameter.

Returns:

The type of the parameter at the given index.

virtual void setParameterType(int index, csm::param::Type pType)

Sets the type of a parameter by its index.

@description Updates the type of the parameter corresponding to the given index.

Parameters:
  • index – The index of the parameter.

  • pType – The new type for the parameter.

virtual double getParameterCovariance(int index1, int index2) const

Retrieves the covariance between two parameters.

@description Provides the covariance between the parameters corresponding to the given indices.

Parameters:
  • index1 – The index of the first parameter.

  • index2 – The index of the second parameter.

Returns:

The covariance between the two parameters.

virtual void setParameterCovariance(int index1, int index2, double covariance)

Sets the covariance between two parameters.

@description Updates the covariance between the parameters corresponding to the given indices.

Parameters:
  • index1 – The index of the first parameter.

  • index2 – The index of the second parameter.

  • covariance – The new covariance value.

virtual int getNumGeometricCorrectionSwitches() const

Retrieves the number of geometric correction switches.

@description Provides the total number of geometric correction switches in the sensor model. SAR sensor models do not use geometric corrections.

Returns:

0, indicating there are no geometric correction switches.

virtual std::string getGeometricCorrectionName(int index) const

Retrieves the name of a geometric correction by its index.

@description Provides the name of the geometric correction corresponding to the given index. SAR sensor models do not use geometric corrections, so this function throws an error if called.

Parameters:

index – The index of the geometric correction.

Returns:

This function will throw an INDEX_OUT_OF_RANGE error.

virtual void setGeometricCorrectionSwitch(int index, bool value, csm::param::Type pType)

Sets a geometric correction switch by its index.

@description Updates the state of the geometric correction switch corresponding to the given index. SAR sensor models do not use geometric corrections, so this function throws an error if called.

Parameters:
  • index – The index of the geometric correction switch.

  • value – The new state for the switch.

  • pType – The parameter type for the switch.

virtual bool getGeometricCorrectionSwitch(int index) const

Retrieves the state of a geometric correction switch by its index.

@description Provides the state of the geometric correction switch corresponding to the given index. SAR sensor models do not use geometric corrections, so this function throws an error if called.

Parameters:

index – The index of the geometric correction switch.

Returns:

This function will throw an INDEX_OUT_OF_RANGE error.

virtual std::vector<double> getCrossCovarianceMatrix(const csm::GeometricModel &comparisonModel, csm::param::Set pSet = csm::param::VALID, const csm::GeometricModel::GeometricModelList &otherModels = csm::GeometricModel::GeometricModelList()) const

Retrieves the cross-covariance matrix between this model and another.

@description This function computes the cross-covariance matrix between the current sensor model and another given model. If the given model is the same as the current model, it returns the covariance matrix of the current model’s parameters. If the models are different, it returns a zero matrix, assuming no correlation between the models.

Parameters:
  • comparisonModel – The geometric model to compare with.

  • pSet – The parameter set for which the covariance is requested.

  • otherModels – List of other geometric models considered in the covariance computation.

Returns:

Vector of doubles representing the flattened cross-covariance matrix.

virtual csm::Version getVersion() const

Returns the version of the sensor model.

@description Provides the version information of the sensor model implementation.

Returns:

Version object containing the major, minor, and revision numbers.

virtual std::string getModelName() const

Gets the model name.

@description Returns the name of the sensor model.

Returns:

String representing the sensor model name.

virtual std::string getPedigree() const

Returns the pedigree of the sensor model.

@description Provides a string that identifies the lineage or source of the sensor model implementation.

Returns:

String representing the sensor model’s pedigree.

virtual std::string getImageIdentifier() const

Retrieves the image identifier.

@description Returns the identifier of the image associated with the sensor model.

Returns:

String representing the image identifier.

virtual void setImageIdentifier(const std::string &imageId, csm::WarningList *warnings = NULL)

Sets the image identifier.

@description Assigns an identifier to the image associated with the sensor model.

Parameters:
  • imageId – The new image identifier.

  • warnings – Optional pointer to a warning list to collect any warnings.

virtual std::string getSensorIdentifier() const

Retrieves the sensor identifier.

@description Returns the identifier of the sensor associated with the sensor model.

Returns:

String representing the sensor identifier.

virtual std::string getPlatformIdentifier() const

Retrieves the platform identifier.

@description Returns the identifier of the platform (satellite/aircraft) associated with the sensor model.

Returns:

String representing the platform identifier.

virtual std::string getCollectionIdentifier() const

Retrieves the collection identifier.

@description Returns a string identifier for the collection to which the image belongs.

Returns:

String representing the collection identifier.

virtual std::string getTrajectoryIdentifier() const

Retrieves the trajectory identifier.

@description Returns a string identifier for the trajectory associated with the sensor model.

Returns:

String representing the trajectory identifier.

virtual std::string getSensorType() const

Gets the sensor type.

@description Returns a string indicating the type of sensor (e.g., SAR for Synthetic Aperture Radar).

Returns:

String representing the sensor type.

virtual std::string getSensorMode() const

Gets the sensor mode.

@description Returns a string indicating the operating mode of the sensor (e.g., “STRIP” for strip mapping).

Returns:

String representing the sensor mode.

virtual std::string getReferenceDateAndTime() const

Gets the reference date and time for the sensor model.

@description Returns a string representing the reference date and time associated with the sensor model.

Returns:

String representing the reference date and time.

virtual csm::Ellipsoid getEllipsoid() const

Retrieves the ellipsoid associated with the sensor model.

@description Returns the ellipsoid used by the sensor model, defined by its semi-major and semi-minor axes.

Returns:

Ellipsoid object representing the ellipsoid used by the sensor model.

virtual void setEllipsoid(const csm::Ellipsoid &ellipsoid)

Sets the ellipsoid for the sensor model.

@description Assigns a new ellipsoid to be used by the sensor model, defined by its semi-major and semi-minor axes.

Parameters:

ellipsoid – The new ellipsoid to use.

void determineSensorCovarianceInImageSpace(csm::EcefCoord &gp, double sensor_cov[4]) const

Determines the sensor covariance in image space.

@description Computes the sensor model covariance matrix transformed into image space for a given ground point.

Parameters:
  • gp – The ground point for which to determine the covariance.

  • sensor_cov – Array to store the resulting 2x2 covariance matrix in image space.

double dopplerShift(csm::EcefCoord groundPt, double tolerance, std::vector<double> adj) const

Calculates the Doppler shift frequency for a given ground point.

@description This function computes the Doppler shift frequency, which is necessary to find the time of the closest approach of the SAR satellite to the given ground point. The Doppler shift is calculated using the relative velocity between the spacecraft and the ground point, the wavelength of the radar signal, and the slant range distance.

Parameters:
  • groundPt – The ground point in ECEF coordinates for which the Doppler shift is calculated.

  • tolerance – The tolerance within which the root-finding algorithm should operate to find the time of closest approach.

  • adj – Adjustments to be applied to the sensor model state.

Returns:

The time at which the Doppler shift frequency is minimized, indicating the closest approach to the ground point.

double slantRange(csm::EcefCoord surfPt, double time, std::vector<double> adj) const

Calculates the slant range distance between the sensor and a ground point.

@description The slant range is the straight-line distance from the satellite to the ground point at a given time.

Parameters:
  • surfPt – The ground point in ECEF coordinates.

  • time – The time at which the slant range is calculated.

  • adj – Adjustments to be applied to the sensor model state.

Returns:

The slant range distance between the satellite and the ground point.

double slantRangeToGroundRange(const csm::EcefCoord &groundPt, double time, double slantRange, double tolerance) const

Converts slant range distance to ground range distance.

@description This function finds the ground range distance corresponding to a given slant range distance for a SAR satellite at a specific time. The ground range is the projection of the slant range onto the ground (Earth’s surface).

Parameters:
  • groundPt – The ground point in ECEF coordinates.

  • time – The time at which the conversion is performed.

  • slantRange – The slant range distance to be converted.

  • tolerance – The tolerance within which the solution for the ground range should be found.

Returns:

The ground range distance corresponding to the given slant range.

double groundRangeToSlantRange(double groundRange, const std::vector<double> &coeffs) const

Converts ground range distance to slant range distance.

@description This function calculates the slant range distance from a given ground range distance using a polynomial defined by the range coefficients. The slant range is the straight-line distance from the satellite to a point on the ground, whereas the ground range is the distance along the ground.

Parameters:
  • groundRange – The ground range distance to be converted to slant range.

  • coeffs – The coefficients of the polynomial used for the conversion.

Returns:

The slant range distance corresponding to the given ground range.

csm::EcefVector getSunPosition(const double imageTime) const

Calculates the sun position at a given image time.

@description Determines the sun’s position in ECEF coordinates at a specified image time. If multiple sun positions are available, Lagrange interpolation is applied. If only one position is available but velocities are provided, linear extrapolation is used. Otherwise, the single provided position is returned.

Parameters:

imageTime – The time of the image for which to calculate the sun position.

Returns:

ECEF vector representing the sun’s position.

std::vector<double> getRangeCoefficients(double time) const

Retrieves range coefficients interpolated for a specific time.

@description Computes the range coefficients at a given time by interpolating the provided scale conversion coefficients. If there are multiple sets of coefficients, Lagrange interpolation is used based on the scale conversion times. Otherwise, the single set of coefficients is returned directly.

Parameters:

time – The time at which to compute the range coefficients.

Returns:

Vector of interpolated range coefficients.

double getValue(int index, const std::vector<double> &adjustments) const

Retrieves the value of a parameter with adjustments applied.

@description Calculates the value of the specified parameter index by adding the corresponding adjustment from the adjustments vector to the parameter’s current value.

Parameters:
  • index – The index of the parameter.

  • adjustments – Vector of adjustments to be applied to the parameters.

Returns:

The adjusted value of the parameter.

Public Members

csm::NoCorrelationModel _NO_CORR_MODEL
std::vector<double> _NO_ADJUSTMENT
std::string m_imageIdentifier
std::string m_platformName
std::string m_sensorName
int m_nLines
int m_nSamples
double m_exposureDuration
double m_scaledPixelWidth
double m_startingEphemerisTime
double m_centerEphemerisTime
double m_endingEphemerisTime
double m_majorAxis
double m_minorAxis
std::string m_platformIdentifier
std::string m_sensorIdentifier
std::string m_trajectoryIdentifier
std::string m_collectionIdentifier
double m_refElevation
double m_minElevation
double m_maxElevation
double m_dtEphem
double m_t0Ephem
std::vector<double> m_scaleConversionCoefficients
std::vector<double> m_scaleConversionTimes
std::vector<double> m_positions
std::vector<double> m_velocities
std::vector<double> m_currentParameterValue
std::vector<csm::param::Type> m_parameterType
csm::EcefCoord m_referencePointXyz
std::vector<double> m_covariance
std::vector<double> m_sunPosition
std::vector<double> m_sunVelocity
double m_wavelength
LookDirection m_lookDirection
std::vector<double> m_noAdjustments
std::shared_ptr<spdlog::logger> m_logger = spdlog::get("usgscsm_logger")

Public Static Functions

static std::string getModelNameFromModelState(const std::string &model_state)

Retrieves the model name from the sensor model state JSON.

@description Parses the model state JSON string to extract the sensor model name. Throws an error if the model name is missing from the state or if the model is not supported.

Parameters:

model_state – The sensor model state in JSON string format.

Throws:

csm::Error – if “m_modelName” is missing or if the sensor model is not supported.

Returns:

The name of the sensor model.

static void applyTransformToState(ale::Rotation const &r, ale::Vec3d const &t, std::string &stateString)

Applies a geometric transformation to the sensor model state.

@description This function applies a rotation and translation to the positions and velocities within the sensor model’s state. The transformation is also applied to the model’s reference point. Note that the Sun’s position and velocity are not altered by this transformation due to their vast distance.

Parameters:
  • r – The rotation to be applied.

  • t – The translation to be applied.

  • stateString – The JSON string representation of the model’s state to which the transformation will be applied. This string is modified in place.

Public Static Attributes

static const std::string _SENSOR_MODEL_NAME = "USGS_ASTRO_SAR_SENSOR_MODEL"
static const int NUM_PARAM_TYPES = 4
static const std::string PARAM_STRING_ALL[] = {"NONE", "FICTITIOUS", "REAL", "FIXED"}
static const csm::param::Type PARAM_CHAR_ALL[] = {csm::param::NONE, csm::param::FICTITIOUS, csm::param::REAL, csm::param::FIXED}
static const int NUM_PARAMETERS = 6
static const std::string PARAMETER_NAME[] = {"X Pos. Bias   ", "Y Pos. Bias   ", "Z Pos. Bias   ", "X Vel. Bias   ", "Y Vel. Bias   ", "Z Vel. Bias   "}