USGSCSM Astro Linescan Sensor Model Class¶
-
class UsgsAstroLsSensorModel : public csm::RasterGM, public virtual csm::SettableEllipsoid¶
- #include <UsgsAstroLsSensorModel.h>
Copyright © 2017-2022 BAE Systems Information and Electronic Systems Integration Inc.
Redistribution and use in source and binary forms, with or without modification, are permitted provided that the following conditions are met:
Redistributions of source code must retain the above copyright notice, this list of conditions and the following disclaimer.
Redistributions in binary form must reproduce the above copyright notice, this list of conditions and the following disclaimer in the documentation and/or other materials provided with the distribution.
Neither the name of the copyright holder nor the names of its contributors may be used to endorse or promote products derived from this software without specific prior written permission.
THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS “AS IS” AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
Public Functions
-
void setState(const std::string &state)¶
-
virtual void replaceModelState(const std::string &stateString)¶
Replaces the sensor model’s state with a new state. @description This function updates the internal state of the sensor model based on a given state string. The state string is parsed, and its contents are used to update model parameters. This function is primarily used for sensor model instantiation or updating the model state in iterative optimization processes.
- Parameters:
stateString – A JSON string representing the new state of the sensor model. This string should contain all necessary model parameters and metadata to fully define the sensor model state.
-
std::string constructStateFromIsd(const std::string imageSupportData, csm::WarningList *list)¶
Constructs the sensor model state from ISD (Image Support Data).
@description Parses the ISD, extracts relevant information, and populates the model state necessary for sensor model operations. This includes sensor, platform, and image identifiers, along with geometric and optical properties. Warnings are collected if any issues arise during parsing.
- Parameters:
imageSupportData – The ISD in string format.
warnings – A pointer to a list for recording any warnings.
- Returns:
A string representation of the sensor model state.
-
void reset()¶
Resets the sensor model to its default state. @description This method reinitializes the sensor model to default values, effectively erasing any adjustments or parameter modifications previously made. It sets all member variables to their default or zero-values, clears any data structures such as vectors, and prepares the model for a fresh start. This includes setting the sensor name to a default value, defining the default focal length, initializing distortion coefficients, and setting various model parameters like the number of lines, samples, and geometric parameters like the major and minor axes of the ellipsoid model used for the planet.
-
UsgsAstroLsSensorModel()¶
Constructor for UsgsAstroLsSensorModel. @description Initializes a new instance of the UsgsAstroLsSensorModel class, setting up the necessary structures and preparing the model for use. This includes allocating resources and setting initial parameter values to their defaults by invoking the reset method.
-
~UsgsAstroLsSensorModel()¶
Destructor for UsgsAstroLsSensorModel. @description Cleans up an instance of UsgsAstroLsSensorModel when it is no longer needed. This includes releasing any resources allocated during the lifetime of the model instance and flushing the logger to ensure all messages are written out. This method ensures a clean shutdown of the sensor model instance, avoiding memory leaks and ensuring that all resources are properly released.
-
virtual std::string getModelState() const¶
Serializes the sensor model state to a JSON string. @description This method captures the current state of the sensor model, including all relevant parameters and configurations, and serializes this information into a structured JSON string.
- Returns:
std::string A JSON string representing the serialized state of the sensor model. The string includes key-value pairs for various model parameters, ensuring a comprehensive snapshot of the model’s current configuration and state.
-
void set(const std::string &state_data)¶
-
virtual csm::ImageCoord groundToImage(const csm::EcefCoord &groundPt, double desiredPrecision = 0.001, double *achievedPrecision = NULL, csm::WarningList *warnings = NULL) const¶
Converts a ground point to its corresponding image coordinates without adjustments.
@description This function estimates the image coordinates (line, sample) for a given ground point in Earth-Centered, Earth-Fixed (ECEF) coordinates. It utilizes the groundToImage function with an internal version that incorporates adjustments. This wrapper function calls the more detailed groundToImage function without providing any adjustments, thereby using the default model parameters.
- Parameters:
ground_pt – The ground point in ECEF coordinates to be converted to image coordinates.
desired_precision – The desired precision for the estimated image coordinates. This influences the iterative process’s termination condition.
achieved_precision – A pointer to a double where the function stores the achieved precision of the estimate. It allows the caller to assess the accuracy of the estimation.
warnings – A pointer to a list for recording any warnings encountered during the computation. It can include warnings about precision not being met.
- Returns:
csm::ImageCoord The estimated image coordinates (line, sample) corresponding to the given ground point.
-
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 its corresponding image coordinates with covariance.
@description This function estimates the image coordinates for a given ground point, including the propagation of the ground point’s covariance through the transformation. It begins by converting the ECEF ground point to image coordinates using the groundToImage method. It then computes the partial derivatives of the line and sample coordinates with respect to the ground point’s X, Y, and Z coordinates. These partial derivatives are used to propagate the ground point’s covariance to the image space, taking into account both modeled and unmodeled errors as well as sensor-specific covariance.
- Parameters:
groundPt – The ground point in ECEF coordinates with covariance.
desired_precision – The desired precision for the image coordinate estimation.
achieved_precision – A pointer to a double where the achieved precision of the estimate is stored.
warnings – A pointer to a list that will be populated with warnings, should they occur.
- Returns:
csm::ImageCoordCovar The image coordinates and their covariance for the given ground point.
-
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 at a given height.
@description This function computes the ground point in Earth-Centered Earth-Fixed (ECEF) coordinates corresponding to a given image coordinate and a specified height above the ellipsoid model of the Earth.
- Parameters:
image_pt – The image coordinate (line and sample) for which to compute the ground point.
height – The height above the ellipsoid at which to compute the ground point, in meters.
desired_precision – The desired precision for the ground point calculation, in meters.
achieved_precision – 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.
- Returns:
csm::EcefCoord The computed 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 associated uncertainties to a ground point in ECEF coordinates.
@description This function computes the ground point in ECEF coordinates corresponding to a given image coordinate with associated uncertainties. It extends the basic image-to-ground transformation by also considering the uncertainties in the image coordinate and the specified height above the ellipsoid. The function calculates the covariance matrix of the resulting ground point by propagating the uncertainties through the transformation process. This includes uncertainties from the sensor model itself (sensor covariance), as well as any unmodeled errors and the provided height variance.
- Parameters:
image_pt – The image coordinate with covariance (line, sample, and their covariance matrix).
height – The height above the ellipsoid at which to compute the ground point, in meters.
heightVariance – The variance of the height above the ellipsoid.
desired_precision – The desired precision for the ground point calculation, in meters.
achieved_precision – 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.
- Returns:
csm::EcefCoordCovar The computed 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 the imaging locus for a given image point proximate to a specified ground point.
@description This function computes the line-of-sight (imaging locus) for a specific image point that is proximate to a given ground point. The computed imaging locus is the line that is closest to the given ground point and lies along the sensor’s line-of-sight.
- Parameters:
image_pt – The image coordinate (line and sample) for which to compute the imaging locus.
ground_pt – The ground point in ECEF coordinates near which the imaging locus is to be computed.
desired_precision – The desired precision for the computation, in meters.
achieved_precision – 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.
- Returns:
csm::EcefLocus The computed imaging locus, including a point on the locus and the unit direction vector of the locus in ECEF coordinates.
-
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 a given image point.
@description This function calculates the remote imaging locus, which is an infinitely extending line from the sensor through the given image point and into space. It is used to model the direction from which the sensor is viewing the Earth’s surface at the specified image coordinate.
- Parameters:
image_pt – The image coordinate (line and sample) for which to compute the remote imaging locus.
desired_precision – The desired precision for the calculation, in meters. This parameter influences the accuracy of the computed direction vector but is not directly used in the current implementation.
achieved_precision – A pointer to a double where the achieved precision of the computation will be stored. This is 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.
- Returns:
csm::EcefLocus The computed remote imaging locus, including a point on the locus (sensor position) and a direction vector indicating the line of sight from the sensor through the image point.
-
virtual csm::ImageCoord getImageStart() const¶
Retrieves the start coordinate of the image.
- Returns:
The start coordinate of the image.
-
virtual csm::ImageVector getImageSize() const¶
Retrieves the size of the image.
- Returns:
The size of the image.
-
virtual std::pair<csm::ImageCoord, csm::ImageCoord> getValidImageRange() const¶
Retrieves the valid range of image coordinates for this sensor model.
- Returns:
A pair of image coordinates specifying the valid range, where the first element is the start coordinate (top-left corner) and the second element is the end coordinate (bottom-right corner).
-
virtual std::pair<double, double> getValidHeightRange() const¶
Retrieves the valid height range for the sensor model.
- Returns:
A pair representing the minimum and maximum valid heights.
-
virtual csm::EcefVector getIlluminationDirection(const csm::EcefCoord &groundPt) const¶
Retrieves the illumination direction for a given ground point.
- Parameters:
groundPt – The ground point for which the illumination direction is queried.
- Returns:
The unit vector representing the direction of illumination.
-
virtual double getImageTime(const csm::ImageCoord &imagePt) const¶
Computes the image time for a given image point.
@description Calculates the time at which a specific image line (and implicitly, the entire image frame) was acquired.
- Parameters:
image_pt – The image point (line, sample) for which to calculate the image acquisition time.
- Returns:
The image acquisition time for the specified image point.
-
virtual csm::EcefCoord getSensorPosition(const csm::ImageCoord &imagePt) const¶
Retrieves the sensor position for a specified image point.
@description This function calculates the sensor’s position in Earth-Centered Earth-Fixed (ECEF) coordinates at the time the specified image point was acquired.
- Parameters:
imagePt – The image coordinates (line, sample) for which to calculate the sensor position.
- Returns:
The sensor’s position in ECEF coordinates.
-
virtual csm::EcefCoord getSensorPosition(double time) const¶
Retrieves the sensor position for a specified time.
@description Calculates the sensor’s position in Earth-Centered Earth-Fixed (ECEF) coordinates at a specific time relative to the start of the image acquisition.
- Parameters:
time – The time at which to calculate the sensor’s position.
- Returns:
The sensor’s position in ECEF coordinates.
-
virtual csm::EcefVector getSensorVelocity(const csm::ImageCoord &imagePt) const¶
Retrieves the sensor velocity for a specified image point.
@description This function calculates the sensor’s velocity in Earth-Centered Earth-Fixed (ECEF) coordinates at the time the specified image point was acquired.
- Parameters:
imagePt – The image coordinates (line, sample) for which to calculate the sensor velocity.
- Returns:
The sensor’s velocity as a vector in ECEF coordinates.
-
virtual csm::EcefVector getSensorVelocity(double time) const¶
Retrieves the sensor velocity for a specified time.
@description Calculates the sensor’s velocity in Earth-Centered Earth-Fixed (ECEF) coordinates at a specific time relative to the start of the image acquisition.
- Parameters:
time – The time at which to calculate the sensor’s velocity.
- Returns:
The sensor’s velocity as a vector in ECEF coordinates.
-
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 partial derivatives of the image coordinates with respect to a specific sensor model parameter.
@description This function calculates how changes in a specific sensor model parameter affect the image coordinates of a ground point.
- Parameters:
index – The index of the sensor model parameter.
ground_pt – The ground point in ECEF coordinates.
desired_precision – The desired computational precision.
achieved_precision – A pointer to a double where the achieved precision will be stored.
warnings – A pointer to a csm::WarningList for logging any warnings.
- Returns:
csm::RasterGM::SensorPartials The partial derivatives of the image line and sample with respect to the specified 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¶
Overloaded method to compute sensor partials using image coordinates.
@description This function calculates how changes in a specific sensor model parameter affect the image coordinates of a ground point; starts with known image coordinates, avoiding the need to recompute them.
- Parameters:
index – The index of the sensor model parameter.
image_pt – The image coordinates.
ground_pt – The ground point in ECEF coordinates.
desired_precision – The desired computational precision.
achieved_precision – A pointer to a double where the achieved precision will be stored.
warnings – A pointer to a csm::WarningList for logging any warnings.
- Returns:
csm::RasterGM::SensorPartials The partial derivatives of the image line and sample with respect to the specified parameter.
-
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 all sensor partial derivatives for a ground point.
@description This function calculates the partial derivatives of the image coordinates with respect to all sensor model parameters for a given ground point.
- Parameters:
ground_pt – The ground point in ECEF coordinates.
pSet – The set of parameters for which partials are to be computed.
desired_precision – The desired computational precision.
achieved_precision – A pointer to a double where the achieved precision will be stored.
warnings – A pointer to a csm::WarningList for logging any warnings.
- Returns:
std::vector<csm::RasterGM::SensorPartials> A vector of sensor partials for all parameters.
-
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 sensor partial derivatives for an image point with respect to all parameters.
@description This function calculates the partial derivatives of the image line and sample coordinates with respect to all adjustable parameters in the sensor model for a given image and ground point.
- Parameters:
image_pt – The image point coordinates (line, sample).
ground_pt – The ground point in ECEF coordinates.
pSet – The set of parameters to consider in the computation.
desired_precision – The desired precision for the computation.
achieved_precision – A pointer to store the achieved precision.
warnings – A pointer to a csm::WarningList for logging any warnings.
- Returns:
A vector of csm::RasterGM::SensorPartials representing the partial derivatives for all parameters.
-
virtual std::vector<double> computeGroundPartials(const csm::EcefCoord &groundPt) const¶
Computes partial derivatives of the image coordinates with respect to the ground point coordinates.
@description This function calculates the partial derivatives of the image line and sample (image coordinates) with respect to the ground point X, Y, and Z coordinates (ECEF).
- Parameters:
ground_pt – The ground point in ECEF coordinates for which to compute partial derivatives.
- Returns:
std::vector<double> A vector containing the six partial derivatives: [dLine/dX, dLine/dY, dLine/dZ, dSample/dX, dSample/dY, dSample/dZ].
-
virtual const csm::CorrelationModel &getCorrelationModel() const¶
Retrieves the correlation model for the sensor.
- Returns:
A reference to the correlation model.
-
virtual std::vector<double> getUnmodeledCrossCovariance(const csm::ImageCoord &pt1, const csm::ImageCoord &pt2) const¶
Retrieves the unmodeled cross covariance between two image points.
- Parameters:
pt1 – The first image point.
pt2 – The second image point.
- Returns:
A vector containing the unmodeled cross covariance values.
-
virtual csm::EcefCoord getReferencePoint() const¶
Retrieves the reference point in ECEF coordinates.
- Returns:
The current reference point.
-
virtual void setReferencePoint(const csm::EcefCoord &groundPt)¶
Sets the reference point in ECEF coordinates.
- Parameters:
ground_pt – The new reference point.
-
virtual int getNumParameters() const¶
Retrieves 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 specific sensor model parameter.
@description Given the index of a sensor parameter, this function returns its name.
- Parameters:
index – The index of the parameter whose name is to be retrieved.
- Returns:
The name of the specified parameter.
-
virtual std::string getParameterUnits(int index) const¶
Retrieves the units of a specific sensor model parameter.
@description Given the index of a sensor parameter, this function returns its units.
- Parameters:
index – The index of the parameter whose units are to be retrieved.
- Returns:
The units of the specified parameter.
Checks if the sensor model has parameters that can be shared across multiple instances.
- Returns:
True if parameters can be shared, false otherwise.
Determines if a specific parameter can be shared across multiple sensor model instances.
- Parameters:
index – The index of the parameter.
- Returns:
True if the parameter can be shared, false otherwise.
-
virtual csm::SharingCriteria getParameterSharingCriteria(int index) const¶
Retrieves the criteria for sharing a specific parameter.
- Parameters:
index – The index of the parameter.
- Returns:
The sharing criteria for the specified parameter.
-
virtual double getParameterValue(int index) const¶
Retrieves the value of a specific sensor model parameter.
@description Given the index of a sensor parameter, this function returns its current value.
- Parameters:
index – The index of the parameter whose value is to be retrieved.
- Returns:
The current value of the specified parameter.
-
virtual void setParameterValue(int index, double value)¶
Sets the value of a specific sensor model parameter.
@description Allows updating the value of a specific sensor parameter identified by its index.
- Parameters:
index – The index of the parameter to update.
value – The new value to set for the parameter.
-
virtual csm::param::Type getParameterType(int index) const¶
Retrieves the type of a specific parameter.
- Parameters:
index – The index of the parameter.
- Returns:
The type of the parameter.
-
virtual void setParameterType(int index, csm::param::Type pType)¶
Sets the type for a specific parameter.
- Parameters:
index – The index of the parameter to set.
pType – The parameter type to set.
-
virtual std::shared_ptr<spdlog::logger> getLogger()¶
Retrieves the logger associated with the sensor model.
@description Accessor for the sensor model’s logging mechanism, allowing for logging of operations and errors.
- Returns:
Shared pointer to the logger.
-
virtual void setLogger(std::string logName)¶
Sets the logger for the sensor model.
@description Associates a named logger with the sensor model, enabling logging of messages under a specific log name.
- Parameters:
logName – The name of the logger.
-
virtual double getParameterCovariance(int index1, int index2) const¶
Retrieves the covariance between two sensor model parameters.
@description Given two parameter indices, this function returns the covariance between them.
- 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 sensor model parameters.
@description Allows for setting the covariance between two specific sensor model parameters.
- Parameters:
index1 – The index of the first parameter.
index2 – The index of the second parameter.
covariance – The new covariance value to set.
-
virtual int getNumGeometricCorrectionSwitches() const¶
Retrieves the number of geometric correction switches available in the model.
- Returns:
The number of geometric correction switches.
-
virtual std::string getGeometricCorrectionName(int index) const¶
Retrieves the name of a specific geometric correction switch.
- Parameters:
index – The index of the correction switch.
- Returns:
The name of the geometric correction switch.
-
virtual void setGeometricCorrectionSwitch(int index, bool value, csm::param::Type pType)¶
Sets the state of a specific geometric correction switch.
- Parameters:
index – The index of the correction switch to set.
value – The new state of the correction switch.
pType – The parameter type associated with the correction switch.
-
virtual bool getGeometricCorrectionSwitch(int index) const¶
Retrieves the state of a specific geometric correction switch.
- Parameters:
index – The index of the correction switch.
- Returns:
The state of the geometric correction switch.
-
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 model.
- Parameters:
comparisonModel – The other geometric model to compare against.
pSet – The parameter set for which the covariance is requested.
otherModels – A list of other geometric models in the scene.
- Returns:
A vector representing the cross covariance matrix.
-
virtual csm::Version getVersion() const¶
Retrieves the version of the sensor model.
- Returns:
A csm::Version object representing the version.
-
virtual std::string getModelName() const¶
Retrieves the name of the sensor model.
- Returns:
The sensor model name.
-
virtual std::string getPedigree() const¶
Retrieves the sensor model’s pedigree.
- Returns:
The pedigree of the sensor model.
-
virtual std::string getImageIdentifier() const¶
Retrieves the image identifier.
- Returns:
The image identifier.
-
virtual void setImageIdentifier(const std::string &imageId, csm::WarningList *warnings = NULL)¶
Sets the image identifier.
- Parameters:
imageId – The new image identifier.
warnings – List to populate with any warnings.
-
virtual std::string getSensorIdentifier() const¶
Retrieves the sensor identifier.
- Returns:
The sensor identifier.
-
virtual std::string getPlatformIdentifier() const¶
Retrieves the platform identifier.
- Returns:
The platform identifier.
-
virtual std::string getCollectionIdentifier() const¶
Retrieves the identifier for the collection to which this sensor model belongs.
- Returns:
A string representing the collection identifier.
-
virtual std::string getTrajectoryIdentifier() const¶
Retrieves the trajectory identifier.
@description Returns a string that identifies the trajectory of the sensor. This could include information like mission name, trajectory number, etc.
- Returns:
A string identifying the trajectory.
-
virtual std::string getSensorType() const¶
Retrieves the sensor type.
- Returns:
A string representing the sensor type.
-
virtual std::string getSensorMode() const¶
Retrieves the sensor mode.
- Returns:
A string representing the sensor mode.
-
virtual std::string getReferenceDateAndTime() const¶
Retrieves the reference date and time for the sensor model.
@description This function calculates and returns the reference date and time corresponding to the sensor model.
- Returns:
A string representing the reference date and time.
-
virtual csm::Ellipsoid getEllipsoid() const¶
Retrieves the ellipsoid used by the sensor model.
- Returns:
A csm::Ellipsoid object representing the ellipsoid.
-
virtual void setEllipsoid(const csm::Ellipsoid &ellipsoid)¶
Sets the ellipsoid used by the sensor model.
- Parameters:
ellipsoid – The new ellipsoid to use.
-
void calculateAttitudeCorrection(const double &time, const std::vector<double> &adj, double attCorr[9]) const¶
Calculates the attitude correction matrix for a given time and set of adjustments.
- Parameters:
time – The time at which to calculate the attitude correction.
adj – The adjustments to apply.
attCorr – The array to store the attitude correction matrix.
-
virtual csm::EcefVector getSunPosition(const double imageTime) const¶
Computes the position of the Sun at a given image time.
@description Determines the Sun’s position in ECEF coordinates at a specified time during the image capture.
- Parameters:
imageTime – The time of interest within the image capture timeline.
- Returns:
The ECEF vector representing the Sun’s position.
Public Members
-
std::string m_imageIdentifier¶
-
std::string m_sensorName¶
-
int m_nLines¶
-
int m_nSamples¶
-
int m_platformFlag¶
-
std::vector<double> m_intTimeLines¶
-
std::vector<double> m_intTimeStartTimes¶
-
std::vector<double> m_intTimes¶
-
double m_startingEphemerisTime¶
-
double m_centerEphemerisTime¶
-
double m_detectorSampleSumming¶
-
double m_detectorLineSumming¶
-
double m_startingDetectorSample¶
-
double m_startingDetectorLine¶
-
int m_ikCode¶
-
double m_focalLength¶
-
double m_zDirection¶
-
DistortionType m_distortionType¶
-
std::vector<double> m_opticalDistCoeffs¶
-
double m_iTransS[3]¶
-
double m_iTransL[3]¶
-
double m_detectorSampleOrigin¶
-
double m_detectorLineOrigin¶
-
double m_mountingMatrix[9]¶
-
double m_majorAxis¶
-
double m_minorAxis¶
-
std::string m_referenceDateAndTime¶
-
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¶
-
double m_dtQuat¶
-
double m_t0Quat¶
-
int m_numPositions¶
-
int m_numQuaternions¶
-
std::vector<double> m_positions¶
-
std::vector<double> m_velocities¶
-
std::vector<double> m_quaternions¶
-
std::vector<double> m_currentParameterValue¶
-
std::vector<csm::param::Type> m_parameterType¶
-
csm::EcefCoord m_referencePointXyz¶
-
double m_gsd¶
-
double m_flyingHeight¶
-
double m_halfSwath¶
-
double m_halfTime¶
-
std::vector<double> m_covariance¶
-
int m_imageFlipFlag¶
-
std::vector<double> m_sunPosition¶
-
std::vector<double> m_sunVelocity¶
-
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 a model state JSON string. @description This function extracts and returns the name of the sensor model from a given state JSON string. It is primarily used to validate or utilize the model name stored within a serialized state representation.
- Parameters:
model_state – A JSON string representing the serialized state of a sensor model. This string must contain a “m_modelName” key with the model’s name as its value.
- Throws:
csm::Error – If the “m_modelName” key is missing from the model state or if the model name in the state is not supported, indicating either an invalid or incompatible state string.
- Returns:
std::string The name of the model extracted from the provided state.
-
static void applyTransformToState(ale::Rotation const &r, ale::Vec3d const &t, std::string &stateString)¶
Applies a rotation and translation transform to the sensor model state. @description This method modifies the sensor model’s state by applying a specified rotation and translation to its elements, including quaternions, positions, and velocities. The operation affects the sensor’s orientation, position, and motion, ensuring that changes in the sensor’s physical setup or perspective are accurately reflected in the model state. The sun’s position and velocity are not modified, assuming their significant distance makes any changes inconsequential to the model.
- Parameters:
r – The rotation to be applied, encapsulated in an
ale::Rotation
object, affecting the sensor’s orientation and angular velocity.t – The translation to be applied, represented as an
ale::Vec3d
vector, shifting the sensor’s position.stateString – A reference to a string holding the serialized JSON representation of the sensor model’s state. This string is modified in-place to reflect the applied transformation.
Public Static Attributes
-
static const std::string _SENSOR_MODEL_NAME = "USGS_ASTRO_LINE_SCANNER_SENSOR_MODEL"¶
-
static const std::string _STATE_KEYWORD[]¶
-
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 = 16¶
-
static const std::string PARAMETER_NAME[] = {"IT Pos. Bias ", "CT Pos. Bias ", "Rad Pos. Bias ", "IT Vel. Bias ", "CT Vel. Bias ", "Rad Vel. Bias ", "Omega Bias ", "Phi Bias ", "Kappa Bias ", "Omega Rate ", "Phi Rate ", "Kappa Rate ", "Omega Accl ", "Phi Accl ", "Kappa Accl ", "Focal Bias "}¶
Private Functions
-
void determineSensorCovarianceInImageSpace(csm::EcefCoord &gp, double sensor_cov[4]) const¶
Calculates sensor covariance in the image coordinate system.
@description This function calculates the covariance matrix of the sensor in the image space for a given ground point. It utilizes the sensor’s partial derivatives with respect to its parameters to transform the parameter covariance into image space, providing a measure of uncertainty in image coordinates that originates from sensor model parameters.
- Parameters:
gp – The ground point in ECEF coordinates for which to determine sensor covariance.
sensor_cov – An array to store the 2x2 covariance matrix. The array will contain the covariance values in row-major order: [line variance, line-sample covariance, sample-line covariance, sample variance].
-
void updateState()¶
Updates the internal state of the sensor model. @description This method computes and updates several key parameters of the sensor model based on the current model state, including the reference point (the ground point corresponding to the center of the image), the ground sample distance (GSD), the flying height of the sensor, the half swath, and the half time duration of the image capture process.
The reference point is calculated by converting the center pixel of the image to a ground point using the imageToGround method. The GSD is then computed as the distance between adjacent pixels in ground space, providing a measure of the spatial resolution of the image. The flying height is determined by calculating the distance from the sensor position to the reference point. Additionally, the method calculates the half swath, which represents half the width of the imaged area on the ground, and the half time duration, which estimates half the time taken to capture the image.
-
double getValue(int index, const std::vector<double> &adjustments) const¶
Retrieves the value of a parameter, adjusted by the specified adjustments.
- Parameters:
index – The index of the parameter.
adjustments – The adjustments to apply to the parameter value.
- Returns:
The adjusted parameter value.
-
virtual csm::ImageCoord groundToImage(const csm::EcefCoord &groundPt, const std::vector<double> &adjustments, double desiredPrecision = 0.001, double *achievedPrecision = NULL, csm::WarningList *warnings = NULL) const¶
Detailed conversion from a ground point to image coordinates, considering adjustments.
@description Performs a refined computation of the image coordinates that view a specified ground point, taking into account model adjustments. This function applies a projective approximation for an initial guess of the image coordinates and then iteratively refines this guess using the secant method to minimize the error in line estimation. The secant method is a numerical technique for finding a root of a function, which in this context is the difference between the estimated and actual lines corresponding to the ground point. The process involves adjusting the initial estimate of the image coordinates based on the given adjustments (if any), and iteratively improving the estimate to achieve or exceed the desired precision level.
- Parameters:
groundPt – The ECEF coordinates of the ground point to project into image space.
adj – Adjustments to apply to the sensor model parameters for this estimation.
desiredPrecision – The precision goal for the iterative estimation process.
achievedPrecision – Pointer to a double where the achieved precision of the estimate is stored.
warnings – Pointer to a list that will be populated with warnings, should they occur.
- Returns:
csm::ImageCoord The refined image coordinates (line, sample) for the ground point.
-
void reconstructSensorDistortion(double &focalX, double &focalY, const double &desiredPrecision) const¶
-
void getQuaternions(const double &time, double quaternion[4]) const¶
Retrieves the quaternions representing the sensor’s orientation at a given time.
- Parameters:
time – The time at which to retrieve the quaternions.
q – The array to store the quaternion values.
-
void losToEcf(const double &line, const double &sample, const std::vector<double> &adj, double &xc, double &yc, double &zc, double &vx, double &vy, double &vz, double &bodyFixedX, double &bodyFixedY, double &bodyFixedZ) const¶
Transforms a line-of-sight vector in image coordinates to the Earth-Centered Fixed (ECF) coordinate system.
@description Computes the ECF coordinates of the sensor position and velocity, as well as the line-of-sight vector in the body-fixed frame. The function accounts for distortion and corrects the line-of-sight vector using the sensor model parameters and adjustments.
- Parameters:
line – The image line coordinate (zero-based).
sample – The image sample coordinate (zero-based, UL pixel center == (0.5, 0.5)).
adj – Adjustments to the sensor model parameters.
xc – Output sensor X coordinate in ECF.
yc – Output sensor Y coordinate in ECF.
zc – Output sensor Z coordinate in ECF.
vx – Output sensor velocity along the X axis.
vy – Output sensor velocity along the Y axis.
vz – Output sensor velocity along the Z axis.
bodyLookX – Output line-of-sight vector X component in body-fixed frame.
bodyLookY – Output line-of-sight vector Y component in body-fixed frame.
bodyLookZ – Output line-of-sight vector Z component in body-fixed frame.
-
void lightAberrationCorr(const double &vx, const double &vy, const double &vz, const double &xl, const double &yl, const double &zl, double &dxl, double &dyl, double &dzl) const¶
Corrects the line-of-sight vector for light aberration.
@description The correction is based on the relative motion between the sensor and the observed surface. This function calculates the correction vector due to light aberration, which is the apparent displacement of the position of a celestial object from its true position (or geometric position), caused by the motion of the observer.
- Parameters:
vx – Sensor velocity along the X axis.
vy – Sensor velocity along the Y axis.
vz – Sensor velocity along the Z axis.
xl – Initial line-of-sight vector X component.
yl – Initial line-of-sight vector Y component.
zl – Initial line-of-sight vector Z component.
dxl – Output corrected line-of-sight vector X component.
dyl – Output corrected line-of-sight vector Y component.
dzl – Output corrected line-of-sight vector Z component.
-
void losEllipsoidIntersect(const double &height, const double &xc, const double &yc, const double &zc, const double &xl, const double &yl, const double &zl, double &x, double &y, double &z, double &achieved_precision, const double &desired_precision, csm::WarningList *warnings = NULL) const¶
Computes the intersection of a line-of-sight vector with an ellipsoid representing the Earth.
@description This function calculates the intersection point of a ray originating from a given camera position (xc, yc, zc) and directed along a vector (xl, yl, zl) with an ellipsoid defined by the sensor model’s major and minor axes. The function also provides the achieved precision of the calculation compared to the desired precision and warns if the line-of-sight does not intersect the ellipsoid.
- Parameters:
height – The height above the ellipsoid at which to calculate the intersection.
xc – The X coordinate of the camera position in ECF.
yc – The Y coordinate of the camera position in ECF.
zc – The Z coordinate of the camera position in ECF.
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 – Output X coordinate of the intersection point in ECF.
y – Output Y coordinate of the intersection point in ECF.
z – Output Z coordinate of the intersection point in ECF.
achieved_precision – The achieved precision of the intersection calculation.
desired_precision – The desired precision for the intersection calculation.
warnings – A list to hold any warnings generated during the calculation.
-
void getAdjSensorPosVel(const double &time, const std::vector<double> &adj, double &xc, double &yc, double &zc, double &vx, double &vy, double &vz, bool calc_vel = true) const¶
Adjusts the sensor position and velocity based on the model parameters and adjustments.
@description This function calculates the adjusted sensor position (xc, yc, zc) and velocity (vx, vy, vz) at a specified time. It considers both nominal sensor positions/velocities and model parameter adjustments. The function can optionally calculate the sensor velocity.
- Parameters:
time – The time at which to calculate the sensor position and velocity.
adj – A vector of adjustments to the sensor model parameters.
xc – Output adjusted sensor X position in ECF.
yc – Output adjusted sensor Y position in ECF.
zc – Output adjusted sensor Z position in ECF.
vx – Output adjusted sensor X velocity.
vy – Output adjusted sensor Y velocity.
vz – Output adjusted sensor Z velocity.
calc_vel – Flag indicating whether to calculate the velocity.
-
std::vector<double> computeDetectorView(const double &time, const csm::EcefCoord &groundPoint, const std::vector<double> &adj) const¶
Computes the detector view for a ground point with adjustments.
@description Calculates focal plane coordinates (x, y) corresponding to a ground point at a specific time, considering model adjustments.
- Parameters:
time – Time at which the observation is made.
groundPoint – The ground point being observed.
adj – Adjustments to the sensor model parameters.
- Returns:
Vector containing the focal plane coordinates (x, y).
-
void computeProjectiveApproximation(const csm::EcefCoord &gp, csm::ImageCoord &ip) const¶
Computes a projective approximation for ground to image conversion.
@description Provides a rapid, approximate conversion from ground coordinates to image coordinates using a projective transformation. If the projective transform has not been initialized, it defaults to the image center.
- Parameters:
gp – The ground point to be approximated.
ip – Reference to store the approximate image point.
-
void createProjectiveApproximation()¶
Creates a projective transformation approximation.
@description Calculates coefficients for a projective transformation based on selected ground to image point mappings, enabling rapid approximate conversions from ground points to image points.
-
double calcDetectorLineErr(double t, csm::ImageCoord const &approxPt, const csm::EcefCoord &groundPt, const std::vector<double> &adj) const¶
Calculates the line error for detector positioning.
@description A helper function for ground to image transformations, which computes the error in detector line positioning relative to an approximated image point and a given ground point.
- Parameters:
t – The time step for adjustment.
approxPt – The approximated image point.
groundPt – The target ground point.
adj – Adjustments to the sensor model parameters.
- Returns:
The calculated line error.