USGSCSM Astro Frame Sensor Model Class

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

Public Functions

UsgsAstroFrameSensorModel()

Initializes the sensor model with default values and logs the creation of the model.

~UsgsAstroFrameSensorModel()
bool isValidModelState(const std::string &stringState, csm::WarningList *warnings)

Checks if the provided model state string represents a valid state for the sensor model. @description Validates the JSON-encoded sensor model state string by checking for the presence of all required keywords and comparing the model name against the expected name.

Parameters:
  • stringState – A JSON-encoded string representing the state of the sensor model.

  • warnings – A pointer to a csm::WarningList for logging any issues found during the validation process. This includes missing required keywords and incorrect model names.

Returns:

A boolean value indicating whether the provided model state string is valid. Returns true if all required keywords are present and the model name matches the expected name; otherwise, returns false.

bool isValidIsd(const std::string &stringIsd, csm::WarningList *warnings)

Validates the Imaging Support Data (ISD) for constructing a sensor model. @description Checks if the provided ISD string can be converted into a valid sensor model state. This serves as a preliminary validation step before attempting to use the ISD for sensor model construction, ensuring that the ISD contains the necessary information.

Parameters:
  • Isd – A JSON-encoded string representing the ISD.

  • warnings – A pointer to a csm::WarningList for logging any issues found during the validation of the ISD.

Returns:

A boolean indicating whether the ISD is valid for model construction. True if a valid model state can be constructed from the ISD, false otherwise.

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

UsgsAstroFrameSensorModel::groundToImage.

Parameters:
  • groundPt

  • desiredPrecision

  • achievedPrecision

  • warnings

Returns:

Returns <line, sample> coordinate in the image corresponding to the ground point without bundle adjustment correction.

std::string constructStateFromIsd(const std::string &jsonIsd, csm::WarningList *warnings)

Constructs the sensor model state from an Imaging Support Data (ISD) string. @description Parses a JSON-encoded ISD string to extract necessary parameters and constructs the state of the sensor model.

Parameters:
  • jsonIsd – A JSON-encoded string representing the ISD from which the state is to be constructed.

  • warnings – A pointer to a csm::WarningList for logging any issues encountered during the parsing of the ISD or the construction of the state.

Throws:

csm::Error – If the ISD does not contain required information or if any errors occur during the parsing process, indicating that the ISD is invalid for creating the sensor model.

Returns:

A JSON-encoded string representing the constructed state of the sensor model.

void reset()

@description Sets the model’s properties and parameters to default values. This includes initializing the sensor’s position, orientation, and other imaging properties, as well as setting up parameter covariances and defaults for imaging geometry.

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

Converts ground coordinates with covariance to image coordinates with covariance. @description Computes the image coordinates corresponding to a given ground point (with covariance), attempting to achieve the specified precision for the conversion. This method is an extension of groundToImage that also considers the covariance of the input ground point but currently returns only the image coordinates without covariance due to partial implementation.

Parameters:
  • groundPt – The ground point (with covariance) in ECEF coordinates to be converted to image coordinates.

  • desiredPrecision – The desired precision for the computation of the image coordinates.

  • achievedPrecision – A pointer to a double where the achieved precision of the computation will be stored. This value indicates how close the computed image point is to the true image point, given the model’s accuracy.

  • warnings – A pointer to a csm::WarningList where any warnings generated during the computation will be added. This can include warnings about the computation process, precision achievements, or any other relevant warnings that do not necessarily indicate a failure.

Returns:

A csm::ImageCoordCovar object representing the computed image coordinates. Due to the partial implementation, the covariance part of the returned object is not populated and should not be considered valid.

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

UsgsAstroFrameSensorModel::groundToImage.

Parameters:
  • groundPt

  • adjustments

  • desired_precision

  • achieved_precision

  • warnings

Returns:

Returns <line,sample> coordinate in the image corresponding to the ground point. This function applies bundle adjustments to the final value.

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

Transforms image coordinates to ground coordinates. @description Converts a point from image coordinates (line, sample) to Earth-Centered, Earth-Fixed (ECEF) ground coordinates at a specified height above the ellipsoid model of the earth. This method accounts for sensor position, orientation, and distortion to accurately map the image point to ground coordinates. Additionally, it can apply corrections for rolling shutter effects.

This function determines if a sample, line intersects the target body and if so, where this intersection occurs in body-fixed coordinates.

Parameters:
  • sample – Sample of the input image.

  • line – Line of the input image.

  • height – ???

  • imagePt – The image coordinates (line and sample) to be converted to ground coordinates.

  • height – The height above the ellipsoid at which to intersect the line of sight from the image point.

  • desiredPrecision – The desired precision for the ground coordinate calculation.

  • achievedPrecision – A pointer to a double where the achieved precision of the calculation will be stored.

  • warnings – A pointer to a csm::WarningList for logging any warnings that occur during the computation.

Returns:

vector<double> Returns the body-fixed X,Y,Z coordinates of the intersection. If no intersection, returns a 3-element vector of 0’s.

Returns:

An ECEF coordinate corresponding to the provided image point.

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

Transforms image coordinates with covariance to ground coordinates with covariance. @description This overload of imageToGround is intended to convert image coordinates with associated covariance to ground coordinates with covariance. Currently, this is an incomplete implementation intended to test software integration capabilities.

Parameters:
  • imagePt – The image coordinates with covariance to be converted to ground coordinates with covariance.

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

  • heightVariance – The variance of the height, contributing to the covariance of the ground point.

  • desiredPrecision – The desired precision for the ground coordinate calculation.

  • achievedPrecision – A pointer to a double where the achieved precision of the calculation will be stored.

  • warnings – A pointer to a csm::WarningList for logging any warnings that occur during the computation.

Returns:

A placeholder ECEF coordinate with covariance indicating an incomplete implementation.

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 the proximate imaging locus for a given image coordinate relative to a specified ground point. @description This function calculates the locus (a line in 3D space) that represents the path of the sensor’s imaging process for a specific image point, adjusted to be as close as possible to a given ground point. The method first computes the remote imaging locus for the image point to determine the imaging direction. It then adjusts the locus to minimize the distance to the specified ground point, effectively creating a proximate imaging locus that represents the imaging geometry relative to that ground point.

Parameters:
  • imagePt – The image coordinate (line and sample) for which to compute the imaging locus.

  • groundPt – The ECEF ground coordinate that the computed locus should be proximate to.

  • desiredPrecision – The desired precision for the computation of the proximate imaging locus.

  • achievedPrecision – A pointer to a double where the achieved precision of the computation will be stored. This value indicates how close the computed locus is to the theoretical perfect locus given the sensor model’s accuracy.

  • warnings – A pointer to a csm::WarningList for logging any warnings that occur during the computation. This can include warnings about the computation process, precision achievements, or any other relevant warnings that do not necessarily indicate a failure.

Returns:

A csm::EcefLocus object representing the computed proximate imaging locus. The returned locus includes both a point on the locus closest to the specified ground point and the direction vector of the locus.

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

Computes the imaging locus for a given image point in remote sensing. @description This function calculates the line of sight (locus) from the sensor model to a point in the image, extending infinitely into space. It determine the direction vector from the sensor’s position through the specifies image point, effectively modeling how the sensor is imaging that point from a distance.

Parameters:
  • imagePt – The image coordinate (line and sample) for which to compute the imaging locus.

  • desiredPrecision – The desired precision for the calculation, which influences the accuracy of the computed direction vector.

  • achievedPrecision – A pointer to a double where the achieved precision of the computation will be stored. This parameter is intended to provide feedback on the computation’s accuracy but is currently set to 0.0 as the calculation is deterministic.

  • warnings – A pointer to a csm::WarningList for logging any warnings that occur during the computation. This might include issues with distortion correction or other computational warnings.

Returns:

A csm::EcefLocus object representing the computed remote imaging locus. The locus includes a point (the sensor’s position) and a direction vector indicating the line of sight from the sensor through the image point.

virtual csm::ImageCoord getImageStart() const

Retrieves the starting coordinates of the image. @description Provides the starting line and sample coordinates of the image as represented by the sensor model. These coordinates indicate the origin (top-left corner) of the image in the sensor’s coordinate system, essential for image processing tasks requiring knowledge of the image’s initial positioning within the sensor’s field of view.

Returns:

csm::ImageCoord containing the starting line and sample of the image. The ‘line’ attribute represents the first line (vertical coordinate), and the ‘samp’ attribute represents the first sample (horizontal coordinate).

virtual csm::ImageVector getImageSize() const

Retrieves the size of the image. @description Returns the total number of lines and samples in the image, providing the image’s dimensions. This information is crucial for various image processing tasks such as cropping, scaling, and geometric corrections, and for allocating appropriate memory for image storage and manipulation.

Returns:

csm::ImageVector containing the size of the image. The ‘line’ attribute represents the total number of lines (vertical extent), and the ‘samp’ attribute represents the total number of samples (horizontal extent).

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

Retrieves the valid image coordinate range of the sensor model. @description Determines the minimum and maximum image coordinates that can be processed by the sensor model. This range is essential for understanding the bounds within which image coordinates are considered valid and can be accurately mapped to ground coordinates.

Returns:

A pair of csm::ImageCoord objects. The first element represents the minimum (starting) image coordinates (line and sample), and the second element represents the maximum image coordinates based on the number of lines and samples in the sensor model.

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

Retrieves the valid height range for the sensor model. @description Provides the minimum and maximum elevation values that can be accurately processed by the sensor model. This range helps in constraining ground point calculations to heights where the model is valid, improving the accuracy of ground-to-image and image-to-ground transformations.

Returns:

A pair of doubles representing the valid height range. The first element is the minimum elevation, and the second element is the maximum elevation, both in meters above the ellipsoid.

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

Calculates the illumination direction for a given ground point. @description Computes the vector representing the direction of illumination at a specific ground point, essential for shading, shadowing, and photometric analysis. The illumination direction is calculated as the vector from the sun to the ground point in the ECEF coordinate system.

Calculates the illumination vector (body-fixed, meters) from the sun to the given ground point.

Parameters:
  • groundPt – The ground point to find the illumination vector to.

  • groundPt – The ECEF coordinates of the ground point for which to compute the illumination direction.

Returns:

csm::EcefVector Returns the illumination vector from the sun to the ground point.

Returns:

A csm::EcefVector representing the direction of illumination at the given ground point. The vector points from the sun towards the ground point, indicating the direction from which the ground point is illuminated.

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

Retrieves the acquisition time of an image point. @description Determines the image acquisition time for a given image coordinate. For this sensor model, the entire image is acquired at once, so the image time for all pixels is the same and set to a constant value.

Parameters:

imagePt – The image coordinate for which to retrieve the acquisition time.

Returns:

The acquisition time of the image point as a double. This model returns 0.0, indicating a single acquisition time for the entire image.

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

Retrieves the sensor’s position at the time of imaging for a given image point. @description Calculates the position of the sensor in ECEF coordinates at the time the specified image point was acquired. For this model, since the entire image is acquired simultaneously, the sensor position is constant for all image points and determined by the current parameter values.

Determines the body-fixed sensor position for the given image coordinate.

Parameters:
  • imagePt – Image coordinate to find the sensor position for.

  • imagePt – The image coordinate for which to determine the sensor’s position.

Throws:
  • csm::Error::BOUNDS – “Image time () out of bounds.”

  • csm::Error – If the given image coordinate is out of the valid range, an error is thrown indicating that the image coordinate is out of bounds.

Returns:

csm::EcefCoord Returns the body-fixed sensor position.

Returns:

The sensor’s position in ECEF coordinates as a csm::EcefCoord object.

virtual csm::EcefCoord getSensorPosition(double time) const

Retrieves the sensor’s position at a specific time. @description Calculates the position of the sensor in ECEF coordinates at the specified time. For this model, the valid image time is a constant (0.0), and the sensor’s position is defined by the current parameter values.

Parameters:

time – The time at which to determine the sensor’s position. For this model, the only valid time is 0.0.

Throws:

csm::Error – If the provided time does not match the valid image time (0.0), an error is thrown indicating that the time is out of bounds.

Returns:

The sensor’s position in ECEF coordinates as a csm::EcefCoord object.

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

Retrieves the sensor’s velocity for a given image point. @description Calculates the velocity of the sensor in ECEF coordinates at the time the specified image point was acquired. For this model, the velocity is constant for all image points and is directly obtained from the ISD.

Determines the velocity of the sensor for the given image coordinate (in body-fixed frame).

Parameters:
  • imagePt – Image coordinate to find the sensor position for.

  • imagePt – The image coordinate for which to determine the sensor’s velocity.

Throws:
  • csm::Error::BOUNDS – “Image coordinate () out of bounds.”

  • csm::Error – If the given image coordinate is out of the valid range, an error is thrown indicating that the image coordinate is out of bounds.

Returns:

csm::EcefVector Returns the sensor velocity in body-fixed frame.

Returns:

The sensor’s velocity in ECEF coordinates as a csm::EcefVector object.

virtual csm::EcefVector getSensorVelocity(double time) const

Retrieves the sensor’s velocity at a specific time. @description Calculates the velocity of the sensor in ECEF coordinates at the specified time. For this model, the only valid image time is 0.0, and the sensor’s velocity is defined by the spacecraft velocity provided in the ISD.

Parameters:

time – The time at which to determine the sensor’s velocity. For this model, the only valid time is 0.0.

Throws:

csm::Error – If the provided time does not match the valid image time (0.0), an error is thrown indicating that the time is out of bounds.

Returns:

The sensor’s velocity in ECEF coordinates as a csm::EcefVector object.

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 with respect to model parameters. @description Calculates the partial derivatives of the image coordinates with respect to a specific model parameter index, given a ground point. This calculation is essential for optimization and analysis tasks where understanding the sensitivity of image coordinates to model parameters is necessary.

Parameters:
  • index – The index of the model parameter for which to compute partials.

  • groundPt – The ground point in ECEF coordinates for which to compute the sensor model partials.

  • desiredPrecision – The desired precision of the partial derivatives calculation.

  • achievedPrecision – A pointer to a double where the achieved precision of the computation will be stored.

  • warnings – A pointer to a csm::WarningList for logging any warnings that occur during the computation of sensor partials.

Returns:

A csm::RasterGM::SensorPartials object representing the partial derivatives of the image coordinate with respect to the specified 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

UsgsAstroFrameSensorModel::computeSensorPartials.

Research: We should investigate using a central difference scheme to approximate the partials. It is more accurate, but it might be costlier calculation-wise.

Parameters:
  • index

  • imagePt

  • groundPt

  • desiredPrecision

  • achievedPrecision

  • warnings

Returns:

The partial derivatives in the line,sample directions.

virtual std::vector<csm::RasterGM::SensorPartials> computeAllSensorPartials(const csm::EcefCoord &groundPt, csm::param::Set pSet = csm::param::VALID, double desiredPrecision = 0.001, double *achievedPrecision = NULL, csm::WarningList *warnings = NULL) const

Computes the sensor model partial derivatives with respect to all model parameters for a ground point. @description Determines the partial derivatives of the image coordinates obtained by projecting a ground point back to the image, with respect to all model parameters. This method first converts the ground point to an image coordinate and then computes the partials for that image point.

Parameters:
  • groundPt – The ground point in ECEF coordinates for which to compute the sensor model partials.

  • pset – The parameter set indicating which model parameters to consider for the partial derivatives calculation.

  • desiredPrecision – The desired precision of the partial derivatives calculation.

  • achievedPrecision – A pointer to a double where the achieved precision of the computation will be stored.

  • warnings – A pointer to a csm::WarningList for logging any warnings that occur during the computation of sensor partials.

Returns:

A vector of csm::RasterGM::SensorPartials objects representing the partial derivatives of the computed image coordinate with respect to each model parameter in the set.

virtual std::vector<csm::RasterGM::SensorPartials> computeAllSensorPartials(const csm::ImageCoord &imagePt, const csm::EcefCoord &groundPt, csm::param::Set pSet = csm::param::VALID, double desiredPrecision = 0.001, double *achievedPrecision = NULL, csm::WarningList *warnings = NULL) const

Computes the sensor model partial derivatives with respect to all model parameters. @description Calculates the partial derivatives of the image coordinates with respect to all model parameters for a given image and ground point pair.

Parameters:
  • imagePt – The image coordinate for which to compute the sensor partials.

  • groundPt – The ground point in ECEF coordinates corresponding to the imagePt.

  • pset – The parameter set indicating which model parameters to consider for the partial derivatives calculation.

  • desiredPrecision – The desired precision of the partial derivatives calculation.

  • achievedPrecision – A pointer to a double where the achieved precision of the computation will be stored.

  • warnings – A pointer to a csm::WarningList for logging any warnings that occur during the computation of sensor partials.

Returns:

A vector of csm::RasterGM::SensorPartials objects representing the partial derivatives of the image coordinate with respect to each model parameter in the set.

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

Computes the partial derivatives of image coordinates with respect to ground coordinates. @description Calculates how changes in the ground point’s ECEF coordinates (X, Y, Z) affect the corresponding image coordinates (line and sample).

Parameters:

groundPt – The ground point in ECEF coordinates for which to compute the partial derivatives of the corresponding image coordinates.

Returns:

A vector of six doubles representing the partial derivatives of image coordinates (line and sample) with respect to ground coordinates (X, Y, Z). The order is [dLine/dX, dLine/dY, dLine/dZ, dSample/dX, dSample/dY, dSample/dZ].

virtual const csm::CorrelationModel &getCorrelationModel() const

Retrieves the correlation model used by the sensor model. @description Returns a reference to the correlation model that describes how image points are correlated with each other in the sensor model. For this implementation, there is no modeled correlation between image points.

Returns:

A constant reference to the sensor model’s CorrelationModel object, which in this case indicates no correlation.

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

Retrieves the unmodeled cross covariance between two image points. @description Computes the cross covariance due to unmodeled errors between two points in the image space. For this implementation, it is assumed that there are no unmodeled errors affecting the covariance between image points.

Parameters:
  • pt1 – The first image coordinate.

  • pt2 – The second image coordinate.

Returns:

A vector of doubles representing the unmodeled cross covariance matrix between the two points, flattened into a vector. Since there are no unmodeled errors, this vector contains only zeros.

virtual csm::Version getVersion() const

Retrieves the version of the sensor model. @description Returns the version of this implementation of the sensor model.

Returns:

A csm::Version object representing the version of the sensor model.

virtual std::string getModelName() const

Retrieves the name of the sensor model. @description Returns the name of the sensor model, which uniquely identifies the model within the CSM framework.

Returns:

A string containing the name of the sensor model.

virtual std::string getPedigree() const

Retrieves the pedigree of the sensor model. @description Returns a string that provides information about the source or lineage of the sensor model, which can be useful for tracking model provenance.

Returns:

A string containing the pedigree of the sensor model.

virtual std::string getImageIdentifier() const

Retrieves the identifier of the image associated with the sensor model. @description Returns the unique identifier of the image that this sensor model represents.

Returns:

A string containing the image identifier.

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

Sets the identifier of the image associated with the sensor model. @description Assigns a new image identifier to the sensor model.

Parameters:
  • imageId – The new image identifier to be set.

  • warnings – A pointer to a csm::WarningList for logging any warnings related to setting the image identifier.

virtual std::string getSensorIdentifier() const

Retrieves the sensor’s identifier. @description Returns the identifier of the sensor that acquired the image.

Returns:

A string containing the sensor’s identifier.

virtual std::string getPlatformIdentifier() const

Retrieves the platform’s identifier. @description Returns the identifier of the platform or spacecraft that hosts the sensor.

Returns:

A string containing the platform’s identifier.

virtual std::string getCollectionIdentifier() const

Retrieves the collection’s identifier. @description Returns the identifier of the image collection to which the image belongs.

Returns:

A string containing the collection’s identifier.

virtual std::string getTrajectoryIdentifier() const

Retrieves the trajectory’s identifier. @description Returns the identifier of the trajectory followed by the platform or sensor. For this model, the trajectory identifier is not used and thus returns an empty string.

Returns:

A string containing the trajectory’s identifier, which is empty for this model.

virtual std::string getSensorType() const

Retrieves the sensor’s type. @description Returns the type of sensor, such as electro-optical, radar, etc., that acquired the image.

Returns:

A string representing the sensor type. For this model, returns CSM_SENSOR_TYPE_EO indicating an electro-optical sensor.

virtual std::string getSensorMode() const

Retrieves the sensor’s mode of operation. @description Returns the mode of operation of the sensor at the time the image was acquired, such as frame, pushbroom, etc.

Returns:

A string representing the sensor mode. For this model, returns CSM_SENSOR_MODE_FRAME, indicating a frame sensor.

virtual std::string getReferenceDateAndTime() const

Retrieves the reference date and time for the sensor model. @description Returns the reference date and time that corresponds to the epoch or start time of the image acquisition.

Returns:

A string representing the reference date and time in a human-readable calendar format.

virtual std::string getModelState() const

Serializes the current state of the sensor model to a string. @description Encodes the complete state of the sensor model into a JSON-formatted string.

Returns:

A string containing the JSON-encoded state of the sensor model. The string starts with the model name for identification, followed by the serialized JSON state with a two-space indentation for readability.

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

Replaces the current model state with a new state. @description Updates the sensor model’s state based on a JSON-encoded state string. This allows for dynamic modification of the sensor model’s parameters and settings, effectively reconfiguring the model.

Parameters:

stringState – A JSON-encoded string representing the new state of the sensor model.

Throws:

csm::Error – If any required state keywords are missing or if the state cannot be successfully parsed and applied, indicating that the new state is invalid or incomplete.

virtual csm::Ellipsoid getEllipsoid() const

Retrieves the ellipsoid associated with the sensor model. @description Provides the ellipsoid used by the sensor model, defined by its major and minor radii.

Returns:

A csm::Ellipsoid object representing the ellipsoid used by the sensor model.

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

Sets the ellipsoid to be used by the sensor model. @description Updates the ellipsoid used by the sensor model to a new one specified by the given csm::Ellipsoid object.

Parameters:

ellipsoid – The new ellipsoid to be used by the sensor model.

virtual csm::EcefCoord getReferencePoint() const

Retrieves the reference point of the sensor model in ECEF coordinates. @description Returns the Earth-Centered, Earth-Fixed (ECEF) coordinates of the reference point.

Returns:

A csm::EcefCoord object representing the reference point in ECEF coordinates.

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

Sets the reference point of the sensor model. @description Updates the sensor model’s reference point to a new location specified by the given ECEF coordinates.

Parameters:

groundPt – The new reference point in ECEF coordinates.

virtual int getNumParameters() const

Retrieves the number of parameters of the sensor model. @description Returns the total number of adjustable parameters within the sensor model.

Returns:

The number of adjustable parameters in the sensor model.

virtual std::string getParameterName(int index) const

Retrieves the name of a specific parameter by its index. @description Provides the name of the parameter at the specified index.

Parameters:

index – The index of the parameter.

Returns:

A string representing the name of the parameter.

virtual std::string getParameterUnits(int index) const

Retrieves the units of a specific parameter by its index. @description Provides the units in which a parameter is measured.

Parameters:

index – The index of the parameter.

Returns:

A string representing the units of the parameter.

virtual bool hasShareableParameters() const

Checks if the sensor model has any parameters that can be shared. @description Determines if any of the sensor model’s parameters are designed to be shareable with other models or entities.

Returns:

A boolean indicating if the model has shareable parameters.

virtual bool isParameterShareable(int index) const

Checks if a specific parameter is shareable. @description Determines if the parameter at the given index is designed to be shareable with other models or entities.

Parameters:

index – The index of the parameter.

Returns:

A boolean indicating if the parameter is shareable.

virtual csm::SharingCriteria getParameterSharingCriteria(int index) const

Retrieves the sharing criteria for a specific parameter. @description Provides the criteria under which a parameter can be shared with other models or entities, if sharing is supported.

Parameters:

index – The index of the parameter.

Returns:

A csm::SharingCriteria object describing the conditions under which the parameter can be shared.

virtual double getParameterValue(int index) const

Retrieves the value of a specific parameter by its index. @description Provides the current value of the parameter identified by the given index.

Parameters:

index – The index of the parameter.

Returns:

The current value of the specified parameter.

virtual void setParameterValue(int index, double value)

Sets the value of a specific parameter by its index. @description Updates the value of the parameter identified by the given index.

Parameters:
  • index – The index of the parameter to be updated.

  • value – The new value for the specified parameter.

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

Retrieves the type of a specific parameter by its index. @description Provides the type classification of the parameter, indicating whether it is adjustable, fixed, or another category defined by the csm::param::Type enumeration.

Parameters:

index – The index of the parameter.

Returns:

The type of the specified parameter as defined by the csm::param::Type enumeration.

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

Sets the type of a specific parameter by its index. @description Updates the type classification of the parameter, allowing for changes in how the parameter is treated by the model or algorithms that use the model.

Parameters:
  • index – The index of the parameter to be updated.

  • pType – The new type for the specified parameter as defined by the csm::param::Type enumeration.

virtual double getParameterCovariance(int index1, int index2) const

Retrieves the covariance between two parameters. @description Provides the covariance value between the parameters identified by the two indexes.

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

  • index2 – The index of the second parameter.

Returns:

The covariance between the specified parameters.

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

Sets the covariance between two parameters. @description Updates the covariance value between the parameters identified by the two indexes.

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

  • index2 – The index of the second parameter.

  • covariance – The new covariance value between the specified parameters.

virtual int getNumGeometricCorrectionSwitches() const

Retrieves the number of geometric correction switches. @description Provides the total number of geometric correction switches available in the sensor model.

Returns:

The number of geometric correction switches in the sensor model.

virtual std::string getGeometricCorrectionName(int index) const

Retrieves the name of a specific geometric correction switch by its index. @description Provides the name of the geometric correction switch identified by the given index. This method is not supported in this model, as geometric corrections are not applicable.

Parameters:

index – The index of the geometric correction switch.

Throws:

csm::Error – If the method is called, indicating that geometric correction switches are not supported by this sensor model.

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

Sets the state of a specific geometric correction switch by its index. @description Updates the state of the geometric correction switch identified by the given index. This method is not supported in this model, as geometric corrections are not applicable.

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

  • value – The new state for the specified geometric correction switch.

  • pType – The parameter type associated with the geometric correction switch.

Throws:

csm::Error – If the method is called, indicating that geometric correction switches are not supported by this sensor model.

virtual bool getGeometricCorrectionSwitch(int index) const

Retrieves the state of a specific geometric correction switch by its index. @description Provides the state of the geometric correction switch identified by the given index. This method is not supported in this model, as geometric corrections are not applicable.

Parameters:

index – The index of the geometric correction switch.

Throws:

csm::Error – If the method is called, indicating that geometric correction switches are not supported by this sensor model.

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

Retrieves the cross-covariance matrix between this model and a comparison model. @description Provides a matrix representing the cross-covariance between the adjustable parameters of this model and those of another geometric model. Since no correlation is assumed between models, a matrix of zeroes is returned.

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

  • pSet – The set of parameters to consider in the covariance computation.

  • otherModels – A list of other geometric models in the scene; not used as no cross-model covariance is supported.

Returns:

A vector of doubles representing the flattened cross-covariance matrix, filled with zeroes.

virtual std::shared_ptr<spdlog::logger> getLogger()

Accessor for the sensor model’s logger. @description Provides access to the sensor model’s logger for logging messages.

Returns:

A shared pointer to the sensor model’s spdlog logger.

virtual void setLogger(std::string logName)

Sets the logger for the sensor model. @description Assigns a logger to the sensor model for logging messages.

Parameters:

logName – The name of the logger to be used by the sensor model.

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

Retrieves a parameter value, optionally adjusted by a given set of adjustments. @description Accesses the value of a parameter identified by its index and applies a given adjustment if provided.

Parameters:
  • index – The index of the parameter whose value is to be retrieved.

  • adjustments – A vector of doubles representing adjustments to the sensor model’s parameters.

Returns:

The adjusted value of the parameter.

void calcRotationMatrix(double m[3][3]) const

Calculates the rotation matrix from the sensor model state or adjusted state. @description Computes the rotation matrix used to rotate vectors from the sensor frame to the earth-centered, earth-fixed (ECEF) frame. This method relies on the quaternion parameters stored in the sensor model’s state.

Parameters:

m – A 3x3 matrix to store the resulting rotation matrix.

void calcRotationMatrix(double m[3][3], const std::vector<double> &adjustments) const

Calculates the rotation matrix given a set of adjustments. @description Computes the rotation matrix using quaternion parameters adjusted by a given set of adjustments.

Parameters:
  • m – A 3x3 matrix to store the resulting rotation matrix.

  • adjustments – A vector of doubles representing the adjustments to the quaternion parameters of the sensor model.

void losEllipsoidIntersect(double height, double xc, double yc, double zc, double xl, double yl, double zl, double &x, double &y, double &z, csm::WarningList *warnings = NULL) const

Computes the intersection of a line-of-sight vector with an ellipsoid. @description This method calculates the point at which a line-of-sight vector, emanating from a given point in space, intersects an ellipsoid model of a planetary body. The ellipsoid is defined by the sensor model’s major and minor axis parameters and can be expanded by a given height to model atmospheric effects or terrain.

Parameters:
  • height – The height above the ellipsoid at which to calculate the intersection.

  • xc – The X-coordinate of the camera position in ECEF space.

  • yc – The Y-coordinate of the camera position in ECEF space.

  • zc – The Z-coordinate of the camera position in ECEF space.

  • xl – The X component of the line-of-sight vector.

  • yl – The Y component of the line-of-sight vector.

  • zl – The Z component of the line-of-sight vector.

  • x – Reference to a double to store the X-coordinate of the intersection point.

  • y – Reference to a double to store the Y-coordinate of the intersection point.

  • z – Reference to a double to store the Z-coordinate of the intersection point.

  • warnings – Optional pointer to a list to store any warnings generated during computation.

Public Members

std::vector<double> m_currentParameterValue
std::vector<double> m_currentParameterCovariance
std::vector<csm::param::Type> m_parameterType
std::vector<double> m_noAdjustments
DistortionType m_distortionType
std::vector<double> m_opticalDistCoeffs
std::vector<double> m_transX
std::vector<double> m_transY
std::vector<double> m_spacecraftVelocity
std::vector<double> m_sunPosition
std::vector<double> m_ccdCenter
std::vector<double> m_iTransS
std::vector<double> m_iTransL
std::vector<double> m_boresight
std::vector<double> m_lineJitter
std::vector<double> m_sampleJitter
std::vector<double> m_lineTimes
double m_majorAxis
double m_minorAxis
double m_focalLength
double m_minElevation
double m_maxElevation
double m_startingDetectorSample
double m_startingDetectorLine
double m_detectorSampleSumming
double m_detectorLineSumming
std::string m_targetName
std::string m_modelName
std::string m_sensorName
std::string m_platformName
std::string m_imageIdentifier
std::string m_collectionIdentifier
double m_ifov
std::string m_instrumentID
double m_focalLengthEpsilon
double m_originalHalfLines
std::string m_spacecraftName
double m_pixelPitch
double m_ephemerisTime
double m_originalHalfSamples
int m_nLines
int m_nSamples
int m_nParameters
csm::EcefCoord m_referencePointXyz

Public Static Functions

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

Extracts the model name from a sensor model’s state string. @description Parses a JSON-encoded string representing the state of a sensor model to extract the model’s name.

Parameters:

model_state – A JSON-encoded string representing the state of the sensor model.

Throws:

csm::Error – If the ‘m_modelName’ key is missing from the model state, indicating an invalid or incomplete state representation, or if the model name extracted does not match the expected sensor model name, indicating that the sensor model is not supported.

Returns:

A string containing the name of the sensor model extracted from the state.

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

Applies a geometric transform to the sensor model’s state. @description Updates the sensor model state by applying a rotation and translation transform.

Parameters:
  • r – The rotation to be applied to the sensor model state, represented as an ale::Rotation object.

  • t – The translation to be applied to the sensor model state, represented as an ale::Vec3d object.

  • stateString – A reference to a string containing the JSON-encoded state of the sensor model. This string will be updated in-place with the transformed state.

Public Static Attributes

static const std::string _SENSOR_MODEL_NAME = "USGS_ASTRO_FRAME_SENSOR_MODEL"
static const int m_numParameters
static const std::string m_parameterName[] = {"X Sensor Position (m)", "Y Sensor Position (m)", "Z Sensor Position (m)", "w", "v1", "v2", "v3"}

Private Members

std::shared_ptr<spdlog::logger> m_logger = spdlog::get("usgscsm_logger")
nlohmann::json _state
csm::NoCorrelationModel _no_corr_model

Private Static Attributes

static const int _NUM_STATE_KEYWORDS
static const int NUM_PARAMETERS = 7
static const std::string _STATE_KEYWORD[]

Friends

friend class UsgsAstroFramePlugin