USGSCSM Utilities Class¶
Functions
-
void computeDistortedFocalPlaneCoordinates(const double &line, const double &sample, const double &sampleOrigin, const double &lineOrigin, const double &sampleSumming, const double &lineSumming, const double &startingSample, const double &startingLine, const double iTransS[], const double iTransL[], double &distortedX, double &distortedY)¶
@description Computes the distorted focal plane coordinates for a given image pixel. This function calculates the position of a pixel in the distorted focal plane given its position in the CCD (Charge-Coupled Device) image sensor coordinate system. The calculation accounts for the CCD coordinate system’s origin, pixel summing, and the starting position of the CCD readout. It uses transformation matrices from the focal plane to CCD samples and lines to perform the conversion. The result is the distorted (x, y) coordinates in the focal plane.
- Parameters:
line – The line (y-coordinate) of the pixel in the CCD coordinate system.
sample – The sample (x-coordinate) of the pixel in the CCD coordinate system.
sampleOrigin – The origin of the CCD coordinate system relative to the top left of the CCD in the x-direction.
lineOrigin – The origin of the CCD coordinate system relative to the top left of the CCD in the y-direction.
sampleSumming – The summing mode of samples (x-direction), indicating how many physical pixels are combined into a single image pixel.
lineSumming – The summing mode of lines (y-direction), indicating how many physical pixels are combined into a single image pixel.
startingSample – The first CCD sample (x-coordinate) corresponding to the start of the image.
startingLine – The first CCD line (y-coordinate) corresponding to the start of the image.
iTransS – An array of three doubles representing the transformation from focal plane to CCD samples. It includes the offset and scaling factors.
iTransL – An array of three doubles representing the transformation from focal plane to CCD lines. It includes the offset and scaling factors.
distortedX – Reference to a double where the x-coordinate of the distorted focal plane position will be stored.
distortedY – Reference to a double where the y-coordinate of the distorted focal plane position will be stored.
-
void removeJitter(const double &line, const double &sample, const std::vector<double> lineJitterCoeffs, const std::vector<double> sampleJitterCoeffs, const std::vector<double> lineTimes, double &dejitteredLine, double &dejitteredSample)¶
@description Computes the de-jittered pixel coordinates given jittered image coordinates, a set of jitter coefficients for both line and sample dimensions, and line exposure times. This function is designed to correct for jitter in images, typically caused by motion or instability during image capture, especially in systems with a rolling shutter. The jitter coefficients should be provided in descending order of power, with the highest power coefficient first and no constant term. For instance, coefficients {1, 2, 3} would correspond to the polynomial 1*t^3 + 2*t^2 + 3*t, where t is the time variable.
- Parameters:
line – The original line (y-coordinate) of the jittered pixel.
sample – The original sample (x-coordinate) of the jittered pixel.
lineJitterCoeffs – A vector of doubles representing the jitter coefficients for the line dimension.
sampleJitterCoeffs – A vector of doubles representing the jitter coefficients for the sample dimension.
lineTimes – A vector of doubles representing the exposure times for each line, used to calculate the actual time at which each line was exposed. This is particularly relevant for rolling shutter cameras where each line of the sensor is exposed at a slightly different time.
dejitteredLine – Reference to a double where the de-jittered line coordinate will be stored.
dejitteredSample – Reference to a double where the de-jittered sample coordinate will be stored.
- Throws:
csm::Error – If the sizes of the jitter coefficient vectors do not match or if the lineTimes vector is empty.
-
void addJitter(const double &line, const double &sample, const double &tolerance, const int &maxIts, const std::vector<double> lineJitterCoeffs, const std::vector<double> sampleJitterCoeffs, const std::vector<double> lineTimes, double &jitteredLine, double &jitteredSample)¶
@description Adds jitter to a given pixel coordinate based on specified jitter coefficients and line exposure times. This function simulates the effect of jitter on image coordinates by applying a reverse process of the jitter removal algorithm. It iteratively applies jitter until the change in the jittered coordinates falls below a specified tolerance or the maximum number of iterations is reached. This approach is useful for testing or simulating the impact of jitter on image data, especially in imaging systems with a rolling shutter.
- Parameters:
line – The original line (y-coordinate) of the pixel before adding jitter.
sample – The original sample (x-coordinate) of the pixel before adding jitter.
tolerance – The tolerance within which the difference between the jittered and original coordinates must fall to consider the jittering process complete.
maxIts – The maximum number of iterations to attempt adding jitter. This prevents infinite loops in cases where the specified tolerance cannot be met.
lineJitterCoeffs – A vector of doubles representing the coefficients for jitter in the line dimension. These coefficients define the polynomial used to model jitter as a function of time.
sampleJitterCoeffs – A vector of doubles representing the coefficients for jitter in the sample dimension. Similar to lineJitterCoeffs, these define the polynomial for modeling jitter.
lineTimes – A vector of doubles representing the exposure times for each line, used to calculate the actual time at which each line was exposed for the purpose of adding jitter based on time.
jitteredLine – Reference to a double where the line coordinate after adding jitter will be stored.
jitteredSample – Reference to a double where the sample coordinate after adding jitter will be stored.
-
void computePixel(const double &distortedX, const double &distortedY, const double &sampleOrigin, const double &lineOrigin, const double &sampleSumming, const double &lineSumming, const double &startingSample, const double &startingLine, const double iTransS[], const double iTransL[], double &line, double &sample)¶
@description Computes the image pixel coordinates from a distorted focal plane coordinate. This function converts coordinates from the distorted focal plane back to the image sensor (CCD) coordinate system, taking into account the CCD’s origin, pixel summing, and starting positions for the sample and line. It uses transformation matrices to map from the focal plane to the CCD coordinates, adjusting for distortion effects that may have been applied previously.
- Parameters:
distortedX – The x-coordinate in the distorted focal plane.
distortedY – The y-coordinate in the distorted focal plane.
sampleOrigin – The x-coordinate of the origin of the CCD coordinate system relative to the top left of the CCD.
lineOrigin – The y-coordinate of the origin of the CCD coordinate system relative to the top left of the CCD.
sampleSumming – The number of physical pixels combined into a single image pixel along the x-direction, affecting resolution and sensitivity.
lineSumming – The number of physical pixels combined into a single image pixel along the y-direction, affecting resolution and sensitivity.
startingSample – The x-coordinate of the first CCD sample corresponding to the start of the image area.
startingLine – The y-coordinate of the first CCD line corresponding to the start of the image area.
iTransS – An array of three doubles representing the transformation from the focal plane to CCD samples, including translation and scaling factors.
iTransL – An array of three doubles representing the transformation from the focal plane to CCD lines, including translation and scaling factors.
line – Reference to a double where the computed line (y-coordinate) in the image pixel coordinate system will be stored.
sample – Reference to a double where the computed sample (x-coordinate) in the image pixel coordinate system will be stored.
-
void calculateRotationMatrixFromQuaternions(double quaternions[4], double cameraToBody[9])¶
@description Calculates a rotation matrix from a quaternion representation. This function converts a given quaternion into its corresponding rotation matrix. Quaternions provide a compact and efficient way to encode the orientation of an object in three-dimensional space. The input quaternion should be normalized or will be normalized within the function to ensure the rotation matrix is orthogonal and properly represents a rotation. The resulting rotation matrix is in row-major order and reflects the rotation in 3D space as described by the quaternion.
- Parameters:
q – An array of four doubles representing the quaternion components in the order: q[0] (w), q[1] (x), q[2] (y), q[3] (z). The quaternion should represent a unit quaternion to correctly encode a rotation.
rotationMatrix – An array of nine doubles where the resulting rotation matrix will be stored. The matrix is stored in row-major order as follows: [ R[0], R[1], R[2], R[3], R[4], R[5], R[6], R[7], R[8] ] where R[] represents the elements of the rotation matrix. This matrix represents the rotation in three-dimensional space according to the quaternion provided.
-
void calculateRotationMatrixFromEuler(double euler[], double rotationMatrix[])¶
@description Calculates a rotation matrix from Euler angles. This function takes in Euler angles (yaw, pitch, and roll) and calculates the corresponding rotation matrix. Euler angles are expected in radians and in the order of rotation about the Z-axis (yaw), Y-axis (pitch), and X-axis (roll). The resulting rotation matrix is in the row-major order and represents the rotation in three-dimensional space according to the right-hand rule.
- Parameters:
euler – An array of three doubles representing the Euler angles in radians:
euler[0]: Rotation angle about the Z-axis (yaw).
euler[1]: Rotation angle about the Y-axis (pitch).
euler[2]: Rotation angle about the X-axis (roll).
rotationMatrix – An array of nine doubles where the resulting rotation matrix will be stored. The matrix is stored in row-major order as follows: [ R[0], R[1], R[2], R[3], R[4], R[5], R[6], R[7], R[8] ] where R[] represents the elements of the rotation matrix.
-
void createCameraLookVector(const double &undistortedFocalPlaneX, const double &undistortedFocalPlaneY, const double &zDirection, const double &focalLength, double cameraLook[])¶
@description Creates a normalized look vector (also known as an imaging ray) in camera space, given coordinates on the undistorted focal plane and the camera’s focal length. This function computes the direction vector from the camera’s optical center through a point on the undistorted focal plane, effectively defining the direction in which the camera is looking at a specific point in the scene. The zDirection parameter allows for adjusting the look vector based on the orientation of the camera’s boresight.
- Parameters:
undistortedFocalPlaneX – The x-coordinate on the undistorted focal plane, typically in units of length (e.g., millimeters), representing the horizontal displacement from the optical axis.
undistortedFocalPlaneY – The y-coordinate on the undistorted focal plane, typically in units of length, representing the vertical displacement from the optical axis.
zDirection – The direction of the camera’s boresight along the z-axis. This should be either 1 or -1, indicating whether the look vector should be oriented along or against the optical axis.
focalLength – The focal length of the camera, defining the distance from the optical center to the focal plane. This value is crucial for determining the depth component of the look vector.
cameraLook – An array of three doubles where the normalized look vector will be stored. The vector is represented as [x, y, z] components in camera space, where the z-axis aligns with the optical axis, and the x and y axes correspond to the horizontal and vertical displacements on the focal plane, respectively.
-
void lagrangeInterp(const int &numTime, const double *valueArray, const double &startTime, const double &delTime, const double &time, const int &vectorLength, const int &i_order, double *valueVector)¶
@description Performs Lagrange interpolation for equally spaced data to estimate the value at a given time point. This function is designed for uniform data intervals and supports up to 8th order interpolation. It gracefully handles points far from the data center to avoid calculation failures. The interpolation order is dynamically adjusted based on the proximity of the desired time to the data points, ensuring optimal accuracy while avoiding overfitting. The function is particularly useful for estimating values from a discrete set of data points in applications such as signal processing or numerical analysis.
- Parameters:
numTime – The number of time points in the dataset.
valueArray – Pointer to the array of values corresponding to each time point. The values for each point are assumed to be stored consecutively in memory.
startTime – The start time of the dataset.
delTime – The uniform time interval between consecutive data points.
time – The time at which the value is to be estimated using interpolation.
vectorLength – The length of the vectors stored in valueArray. This parameter allows the function to handle multidimensional data by treating each set of vectorLength consecutive elements in valueArray as a vector associated with a single time point.
i_order – The maximum order of interpolation desired. The actual order used may be less than this value depending on the proximity of the interpolation point to the data boundaries.
valueVector – Pointer to the array where the interpolated value(s) will be stored. This array must be pre-allocated with at least vectorLength elements. The interpolated values are stored as a vector of length vectorLength.
- Throws:
csm::Error – If numTime is less than 2, indicating insufficient data points for interpolation.
-
double brentRoot(double lowerBound, double upperBound, std::function<double(double)> func, double epsilon = 1e-10)¶
@description Finds a root of a function within a specified interval using Brent’s method. Brent’s method is an efficient algorithm for finding roots of a function. It combines the bisection method, the secant method, and inverse quadratic interpolation. It has the reliability of bisection but can converge much faster when close to the root. The function to find the root of is provided as a parameter, allowing for flexible usage with any continuous function.
- Parameters:
lowerBound – The lower bound of the interval in which to search for the root.
upperBound – The upper bound of the interval in which to search for the root. It is assumed that the function crosses the x-axis at least once in the interval [lowerBound, upperBound].
func – The function for which to find the root. This is a std::function object that takes a double and returns a double. The function should be continuous, and func(lowerBound) and func(upperBound) should have opposite signs.
epsilon – The tolerance for convergence. The algorithm stops when the interval between the current best guess of the root and the previous best guess is less than epsilon, or when the function value at the current guess is within epsilon of zero.
- Throws:
std::invalid_argument – If func(lowerBound) and func(upperBound) have the same sign, indicating that Brent’s method may not be applicable as there might not be a root within the interval or the function is not continuous.
- Returns:
The estimated root of the function within the specified interval and tolerance.
-
void newtonRaphson(double dx, double dy, double &ux, double &uy, std::vector<double> const &opticalDistCoeffs, DistortionType distortionType, const double tolerance, std::function<void(double, double, double&, double&, std::vector<double> const&)> distortionFunction, std::function<void(double, double, double*, std::vector<double> const&)> distortionJacobian)¶
@description Applies the Newton-Raphson method to un-distort a pixel coordinate (dx, dy), producing the undistorted coordinate (ux, uy).
- Parameters:
dx – The distorted x-coordinate of the pixel.
dy – The distorted y-coordinate of the pixel.
ux – Reference to a double where the undistorted x-coordinate will be stored.
uy – Reference to a double where the undistorted y-coordinate will be stored.
opticalDistCoeffs – A constant reference to a vector of doubles containing the optical distortion coefficients. These coefficients are specific to the distortion model being used.
distortionType – An enumeration of the type of distortion model to be corrected. This parameter dictates how the distortion and its Jacobian functions interpret the distortion coefficients.
tolerance – The tolerance level for the convergence of the Newton-Raphson method. The iterative process will stop when the sum of the absolute differences between the distorted coordinates and their corrections falls below this tolerance or when the maximum number of iterations is reached.
distortionFunction – A std::function object representing the distortion function. It takes distorted coordinates (x, y), a reference to the distortion coefficients, and outputs the distortion effects (fx, fy) to be applied to (x, y).
distortionJacobian – A std::function object representing the Jacobian of the distortion function. It takes distorted coordinates (x, y) and outputs the Jacobian matrix as an array of 4 doubles, which is used to calculate the next iteration’s correction.
-
double evaluatePolynomial(const std::vector<double> &coeffs, double x)¶
-
double evaluatePolynomialDerivative(const std::vector<double> &coeffs, double x)¶
@description Evaluates a polynomial at a given point x. The polynomial is defined by its coefficients, with the coefficients provided in ascending order of power. This function calculates the value of the polynomial for a specific x-value using Horner’s method.
- Parameters:
coeffs – A constant reference to a vector of doubles representing the coefficients of the polynomial. The first element represents the coefficient of the lowest power term, and the last element represents the coefficient of the highest power term.
x – The point at which to evaluate the polynomial.
- Throws:
std::invalid_argument – If the coeffs vector is empty, as a polynomial requires at least one coefficient.
- Returns:
The value of the polynomial at the given point x.
-
double polynomialRoot(const std::vector<double> &coeffs, double guess, double threshold = 1e-10, int maxIterations = 30)¶
@description Finds an approximate root of a polynomial near a given initial guess using Newton’s method. This iterative algorithm refines the guess based on the polynomial and its derivative, converging to a root when the value of the polynomial at the current guess is below a specified threshold. Newton’s method is efficient and effective for well-behaved functions and starting guesses close to the actual root. However, convergence is not guaranteed if the function’s derivative is very small or if the initial guess is far from any root.
- Parameters:
coeffs – A vector of doubles representing the coefficients of the polynomial to find the root of. The coefficients should be ordered from the lowest degree term to the highest.
guess – An initial guess for the root of the polynomial.
threshold – The value below which the absolute value of the polynomial is considered close enough to zero to conclude that the root has been found.
maxIterations – The maximum number of iterations to attempt before giving up on finding a root.
- Throws:
std::invalid_argument – If the derivative of the polynomial at the current guess is too close to zero, making it impossible to proceed with the next iteration. This condition indicates either a poor initial guess or that the guess is near a point where the polynomial’s slope is horizontal.
std::invalid_argument – If the function does not converge to a root within the specified number of iterations. This may indicate that the initial guess is too far from any root or that the specified threshold and maximum number of iterations are too strict.
- Returns:
The approximate root of the polynomial, found near the initial guess.
-
double computeEllipsoidElevation(double x, double y, double z, double semiMajor, double semiMinor, double desired_precision = 0.001, double *achieved_precision = nullptr)¶
@description Computes the elevation of a point above or below an ellipsoidal model of the Earth, given the point’s Cartesian coordinates. The function iteratively adjusts the elevation estimate until the change between iterations is less than the specified desired precision or until a maximum number of iterations is reached.
- Parameters:
x – The x-coordinate of the point in Cartesian coordinate system.
y – The y-coordinate of the point in Cartesian coordinate system.
z – The z-coordinate of the point in Cartesian coordinate system.
semiMajor – The semi-major axis of the ellipsoid, typically representing the equatorial radius.
semiMinor – The semi-minor axis of the ellipsoid, typically representing the polar radius.
desired_precision – The desired precision for the elevation calculation. The iterative process will stop once the change in elevation estimate between iterations is less than this value.
achieved_precision – Pointer to a double where the achieved precision of the calculation will be stored. This is the absolute difference in elevation estimate between the final two iterations. If null, this output is not provided.
- Returns:
The elevation of the point above (positive value) or below (negative value) the ellipsoidal model, in the same units as the semi-major and semi-minor axes.
-
csm::EcefVector operator*(double scalar, const csm::EcefVector &vec)¶
@description Multiplies a scalar by an ECEF vector.
- Parameters:
scalar – A double value to scale the ECEF vector.
vec – A constant reference to an ECEF vector to be scaled.
- Returns:
An ECEF vector where each component is multiplied by the scalar.
-
csm::EcefVector operator*(const csm::EcefVector &vec, double scalar)¶
@description Multiplies an ECEF vector by a scalar.
- Parameters:
vec – A constant reference to an ECEF vector to be scaled.
scalar – A double value to scale the ECEF vector.
- Returns:
An ECEF vector where each component is multiplied by the scalar.
-
csm::EcefVector operator/(const csm::EcefVector &vec, double scalar)¶
@description Divides an ECEF vector by a scalar.
- Parameters:
vec – A constant reference to an ECEF vector to be divided.
scalar – A double value by which the ECEF vector components are divided.
- Returns:
An ECEF vector where each component is divided by the scalar.
-
csm::EcefVector operator+(const csm::EcefVector &vec1, const csm::EcefVector &vec2)¶
@description Adds two ECEF vectors.
- Parameters:
vec1 – A constant reference to the first ECEF vector.
vec2 – A constant reference to the second ECEF vector.
- Returns:
An ECEF vector resulting from component-wise addition of vec1 and vec2.
-
csm::EcefVector operator-(const csm::EcefVector &vec1, const csm::EcefVector &vec2)¶
@description Subtracts one ECEF vector from another.
- Parameters:
vec1 – A reference to the ECEF vector from which vec2 is subtracted.
vec2 – A reference to the ECEF vector to subtract from vec1.
- Returns:
An ECEF vector resulting from component-wise subtraction of vec2 from vec1.
-
double dot(const csm::EcefVector &vec1, const csm::EcefVector &vec2)¶
@description Calculates the dot product of two ECEF vectors.
- Parameters:
vec1 – A constant reference to the first ECEF vector.
vec2 – A constant reference to the second ECEF vector.
- Returns:
A double representing the dot product of vec1 and vec2.
-
csm::EcefVector cross(const csm::EcefVector &vec1, const csm::EcefVector &vec2)¶
@description Calculates the cross product of two ECEF vectors.
- Parameters:
vec1 – A constant reference to the first ECEF vector.
vec2 – A constant reference to the second ECEF vector.
- Returns:
An ECEF vector representing the cross product of vec1 and vec2. The resulting vector is orthogonal to both vec1 and vec2.
-
double norm(const csm::EcefVector &vec)¶
@description Calculates the norm (magnitude) of an ECEF vector.
- Parameters:
vec – A constant reference to the ECEF vector.
- Returns:
A double value representing the Euclidean norm of the vector, which is the square root of the sum of the squares of its components.
-
csm::EcefVector normalized(const csm::EcefVector &vec)¶
@description Normalizes an ECEF vector, scaling it to a unit vector.
- Parameters:
vec – A constant reference to the ECEF vector to be normalized.
- Returns:
An ECEF vector that points in the same direction as the original vector but with a norm (magnitude) of 1. If the input vector is a zero vector, the behavior is undefined.
-
csm::EcefVector projection(const csm::EcefVector &vec1, const csm::EcefVector &vec2)¶
@description Calculates the projection of one ECEF vector onto another.
- Parameters:
vec1 – A constant reference to the ECEF vector being projected.
vec2 – A constant reference to the ECEF vector upon which vec1 is projected.
- Returns:
An ECEF vector representing the projection of vec1 onto vec2. This is the component of vec1 that lies in the direction of vec2.
-
csm::EcefVector rejection(const csm::EcefVector &vec1, const csm::EcefVector &vec2)¶
@description Calculates the rejection of one ECEF vector from another.
- Parameters:
vec1 – A constant reference to the ECEF vector being rejected.
vec2 – A constant reference to the ECEF vector from which vec1 is rejected.
- Returns:
An ECEF vector representing the rejection of vec1 from vec2. This is the component of vec1 that is orthogonal to vec2, effectively vec1 minus its projection on vec2.
-
double metric_conversion(double val, std::string from, std::string to = "m")¶
@description Converts a measurement value from one metric unit to another.
- Parameters:
val – The numerical value to be converted.
from – The unit of measure of the input value (e.g., “m” for meters, “km” for kilometers).
to – The target unit of measure for the conversion (e.g., “m” for meters, “km” for kilometers).
- Returns:
The converted value in the target unit of measure.
-
std::string getSensorModelName(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Extracts the sensor model name from an ISD JSON object.
- Parameters:
isd – The ISD JSON object containing sensor model information.
list – Pointer to a csm::WarningList where warnings can be added if the sensor model name cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The sensor model name as a string. Returns an empty string if the name cannot be found or parsed from the ISD.
-
std::string getImageId(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Extracts the image identifier from an ISD JSON object.
- Parameters:
isd – The ISD JSON object containing image identification information.
list – Pointer to a csm::WarningList where warnings can be added if the image identifier cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The image identifier as a string. Returns an empty string if the identifier cannot be found or parsed from the ISD.
-
std::string getSensorName(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Extracts the sensor name from an ISD JSON object.
- Parameters:
isd – The ISD JSON object containing the sensor name information.
list – Pointer to a csm::WarningList where warnings can be added if the sensor name cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The sensor name as a string. Returns an empty string if the name cannot be found or parsed from the ISD.
-
std::string getPlatformName(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Extracts the platform name from an ISD JSON object.
- Parameters:
isd – The ISD JSON object containing the platform name information.
list – Pointer to a csm::WarningList where warnings can be added if the platform name cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The platform name as a string. Returns an empty string if the name cannot be found or parsed from the ISD.
-
std::string getLogFile(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Extracts the log file name from an ISD JSON object.
- Parameters:
isd – The ISD JSON object containing the log file name information.
list – Pointer to a csm::WarningList where warnings can be added if the log file name cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The log file name as a string. Returns an empty string if the log file name cannot be found or parsed from the ISD.
-
int getTotalLines(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the total number of lines in an image from an ISD JSON object.
- Parameters:
isd – The ISD JSON object containing metadata about the image.
list – Pointer to a csm::WarningList where warnings can be added if the total number of lines cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The total number of lines in the image as an integer. Returns 0 if the information cannot be found or parsed from the ISD.
-
int getTotalSamples(nlohmann::json isd, csm::WarningList *list = nullptr)¶
Retrieves the total number of samples in an image from an ISD JSON object.
- Parameters:
isd – The ISD JSON object containing metadata about the image.
list – Pointer to a csm::WarningList where warnings can be added if the total number of samples cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The total number of samples in the image as an integer. Returns 0 if the information cannot be found or parsed from the ISD.
-
double getStartingTime(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the starting ephemeris time of the image capture from an ISD JSON object.
- Parameters:
isd – The ISD JSON object containing temporal data about the image.
list – Pointer to a csm::WarningList where warnings can be added if the starting ephemeris time of the image capture cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The starting ephemeris time of the image capture as a double. Returns 0.0 if the information cannot be found or parsed from the ISD.
-
double getCenterTime(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the center ephemeris time of the image capture from an ISD JSON object.
- Parameters:
isd – The ISD JSON object containing temporal data about the image.
list – Pointer to a csm::WarningList where warnings can be added if the center ephemeris time cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The center ephemeris time of the image capture as a double. Returns 0.0 if the information cannot be found or parsed from the ISD.
-
double getEndingTime(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the ending ephemeris time of the image capture from an ISD JSON object.
- Parameters:
isd – The ISD JSON object containing temporal data about the image.
list – Pointer to a csm::WarningList where warnings can be added if the ending ephemeris time cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The ending ephemeris time of the image capture as a double. Returns 0.0 if the information cannot be found or parsed from the ISD.
-
std::vector<double> getIntegrationStartLines(std::vector<std::vector<double>> lineScanRate, csm::WarningList *list = nullptr)¶
@description Extracts the integration start lines from a line scan rate table.
- Parameters:
lineScanRate – A vector of vectors, each inner vector containing three doubles: the start line, the start time, and the duration of each integration period.
list – Pointer to a csm::WarningList where warnings are added if any issues occur during parsing.
- Returns:
A vector of doubles, each representing the start line of an integration period extracted from the line scan rate table.
-
std::vector<double> getIntegrationStartTimes(std::vector<std::vector<double>> lineScanRate, csm::WarningList *list = nullptr)¶
@description Extracts the integration start times from a line scan rate table.
- Parameters:
lineScanRate – A vector of vectors, each inner vector containing three doubles: the start line, the start time, and the duration of each integration period.
list – Pointer to a csm::WarningList where warnings are added if any issues occur during parsing.
- Returns:
A vector of doubles, each representing the start time of an integration period extracted from the line scan rate table.
-
std::vector<double> getIntegrationTimes(std::vector<std::vector<double>> lineScanRate, csm::WarningList *list = nullptr)¶
@description Extracts the integration times (durations) from a line scan rate table.
- Parameters:
lineScanRate – A vector of vectors, each inner vector containing three doubles: the start line, the start time, and the duration of each integration period.
list – Pointer to a csm::WarningList where warnings are added if any issues occur during parsing.
- Returns:
A vector of doubles, each representing the duration of an integration period extracted from the line scan rate table.
-
double getExposureDuration(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the line exposure duration from an ISD JSON object. The exposure duration is the time, typically in seconds, for which each line of the sensor was exposed during image capture.
- Parameters:
isd – The ISD JSON object containing the line exposure duration.
list – Pointer to a csm::WarningList where warnings are added if the line exposure duration cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The line exposure duration as a double. If the duration cannot be parsed, 0.0 is returned.
-
double getScaledPixelWidth(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the scaled pixel width from an ISD JSON object. The scaled pixel width is the physical width of a pixel on the sensor, typically in meters, after applying any necessary scale factors.
- Parameters:
isd – The ISD JSON object containing the scaled pixel width.
list – Pointer to a csm::WarningList where warnings are added if the scaled pixel width cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The scaled pixel width as a double. If the width cannot be parsed, 0.0 is returned.
-
std::string getLookDirection(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the look direction from an ISD JSON object. The look direction is a string that indicates the general direction in which the sensor was pointing when the image was captured, such as “NADIR”, “FORWARD”, or “BACKWARD”.
- Parameters:
isd – The ISD JSON object containing the look direction.
list – Pointer to a csm::WarningList where warnings are added if the look direction cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The look direction as a std::string. If the look direction cannot be parsed, an empty string is returned.
-
std::vector<double> getScaleConversionCoefficients(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the range conversion coefficients from an ISD JSON object. These coefficients are used in the conversion of instrument data to physical units.
- Parameters:
isd – The ISD JSON object containing the range conversion coefficients.
list – Pointer to a csm::WarningList where warnings are added if the range conversion coefficients cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
A std::vector<double> containing the range conversion coefficients. If the coefficients cannot be parsed, an empty vector is returned.
-
std::vector<double> getScaleConversionTimes(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the range conversion times from an ISD JSON object. These times are used in the conversion of instrument data to physical units and are typically in seconds.
- Parameters:
isd – The ISD JSON object containing the range conversion times.
list – Pointer to a csm::WarningList where warnings are added if the range conversion times cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
A std::vector<double> containing the range conversion times. If the times cannot be parsed, an empty vector is returned.
-
int getSampleSumming(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the sample summing value from an ISD JSON object. This value represents the pixel summing in the sample (horizontal) direction, which is used to reduce the image resolution for various operational reasons.
- Parameters:
isd – The ISD JSON object containing the imaging data.
list – Pointer to a csm::WarningList where warnings are added if the sample summing value cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The sample summing value as an integer. Returns 0 if the information cannot be found or parsed from the ISD.
-
int getLineSumming(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the detector line summing value from an ISD JSON object.
- Parameters:
isd – The ISD JSON object containing metadata about the imaging sensor.
list – Pointer to a csm::WarningList where warnings are added if the line summing value cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The line summing value as an integer, representing the number of detector pixels summed in the line (vertical) direction to form the image. Returns 0 if the information cannot be found or parsed from the ISD.
-
double getFocalLength(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the focal length of the imaging system from an ISD JSON object.
- Parameters:
isd – The ISD JSON object containing the focal length model of the imaging sensor.
list – Pointer to a csm::WarningList where warnings are added if the focal length cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The focal length as a double, usually in millimeters. Returns 0.0 if the focal length cannot be found or parsed from the ISD.
-
double getFocalLengthEpsilon(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the uncertainty of the focal length (focal epsilon) from an ISD JSON object.
- Parameters:
isd – The ISD JSON object containing the focal length model and its uncertainty.
list – Pointer to a csm::WarningList where warnings are added if the focal length uncertainty cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The focal length uncertainty (epsilon) as a double, usually in millimeters. Returns 0.0 if the uncertainty cannot be found or parsed from the ISD.
-
std::vector<double> getFocal2PixelLines(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the transformation from focal plane coordinates to detector lines as a vector of doubles from an ISD JSON object.
- Parameters:
isd – The ISD JSON object containing the transformation data.
list – Pointer to a csm::WarningList where warnings are added if the transformation cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
A std::vector<double> containing the transformation from focal plane coordinates to detector lines. Returns an empty vector if the transformation cannot be found or parsed from the ISD.
-
std::vector<double> getFocal2PixelSamples(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the transformation from focal plane coordinates to detector samples as a vector of doubles from an ISD JSON object.
- Parameters:
isd – The ISD JSON object containing the transformation data.
list – Pointer to a csm::WarningList where warnings are added if the transformation cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
A std::vector<double> containing the transformation from focal plane coordinates to detector samples. Returns an empty vector if the transformation cannot be found or parsed from the ISD.
-
double getDetectorCenterLine(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the detector center line from an ISD JSON object.
- Parameters:
isd – The ISD JSON object containing the detector center information.
list – Pointer to a csm::WarningList where warnings are added if the detector center line cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The detector center line as a double. Returns 0.0 if the center line cannot be found or parsed from the ISD.
-
double getDetectorCenterSample(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the detector center sample from an ISD JSON object.
- Parameters:
isd – The ISD JSON object containing the detector center information.
list – Pointer to a csm::WarningList where warnings are added if the detector center sample cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The detector center sample as a double. Returns 0.0 if the center sample cannot be found or parsed from the ISD.
-
double getDetectorStartingLine(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the starting line of the detector from an ISD JSON object.
- Parameters:
isd – The ISD JSON object containing the starting line information.
list – Pointer to a csm::WarningList where warnings are added if the starting line cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The starting line of the detector as a double. Returns 0.0 if the starting line cannot be found or parsed from the ISD.
-
double getDetectorStartingSample(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the starting sample of the detector from an ISD JSON object.
- Parameters:
isd – The ISD JSON object containing the starting sample information.
list – Pointer to a csm::WarningList where warnings are added if the starting sample cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The starting sample of the detector as a double. Returns 0.0 if the starting sample cannot be found or parsed from the ISD.
-
double getMinHeight(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the minimum height above the reference ellipsoid from an ISD JSON object. The height is converted to meters if necessary.
- Parameters:
isd – The ISD JSON object containing the minimum height and its unit.
list – Pointer to a csm::WarningList where warnings are added if the minimum height cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The minimum height above the reference ellipsoid in meters as a double. Returns 0.0 if the minimum height cannot be found or parsed from the ISD.
-
double getMaxHeight(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the maximum height above the reference ellipsoid from an ISD JSON object. The height is converted to meters if necessary.
- Parameters:
isd – The ISD JSON object containing the maximum height and its unit.
list – Pointer to a csm::WarningList where warnings are added if the maximum height cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The maximum height above the reference ellipsoid in meters as a double. Returns 0.0 if the maximum height cannot be found or parsed from the ISD.
-
double getSemiMajorRadius(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the semimajor radius of the reference ellipsoid from an ISD JSON object. The radius is converted to meters if necessary.
- Parameters:
isd – The ISD JSON object containing the semimajor radius and its unit.
list – Pointer to a csm::WarningList where warnings are added if the semimajor radius cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The semimajor radius of the reference ellipsoid in meters as a double. Returns 0.0 if the semimajor radius cannot be found or parsed from the ISD.
-
double getSemiMinorRadius(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the semiminor radius of the reference ellipsoid from an ISD JSON object. The radius is converted to meters if necessary.
- Parameters:
isd – The ISD JSON object containing the semiminor radius and its unit.
list – Pointer to a csm::WarningList where warnings are added if the semiminor radius cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The semiminor radius of the reference ellipsoid in meters as a double. Returns 0.0 if the semiminor radius cannot be found or parsed from the ISD.
-
DistortionType getDistortionModel(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Converts the distortion model name from the ISD JSON object to the enumeration type. This function searches the “optical_distortion” field for a key that matches a known distortion model name and converts it to the corresponding DistortionType enum. Defaults to DistortionType::TRANSVERSE if no match is found or the model cannot be parsed.
- Parameters:
isd – The ISD JSON object containing the optical distortion information.
list – Pointer to a csm::WarningList where warnings are added if the distortion model name cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The DistortionType enum corresponding to the distortion model found in the ISD. Returns DistortionType::TRANSVERSE by default if the model name cannot be parsed or is not recognized.
-
DistortionType getDistortionModel(int aleDistortionModel, csm::WarningList *list = nullptr)¶
@description Converts an integer representing a distortion model from ALE (Abstraction Library for Ephemerides) to the corresponding DistortionType enum. This function directly maps integers to the DistortionType enumeration based on the ALE distortion model definitions. Defaults to DistortionType::TRANSVERSE if the integer does not match any known model or the model cannot be parsed.
- Parameters:
aleDistortionModel – An integer representing the ALE distortion model.
list – Pointer to a csm::WarningList where warnings are added if the distortion model integer cannot be mapped to a DistortionType enum. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The DistortionType enum corresponding to the ALE distortion model integer. Returns DistortionType::TRANSVERSE by default if the model integer cannot be mapped or is not recognized.
-
std::vector<double> getDistortionCoeffs(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the distortion coefficients for the specified distortion model from an ISD JSON object. This function first determines the distortion model using getDistortionModel() and then extracts the coefficients specific to that model. If the coefficients cannot be parsed, a default set of coefficients is returned based on the model, and a warning is added to the list.
- Parameters:
isd – The ISD JSON object containing the optical distortion coefficients for the specific distortion model.
list – Pointer to a csm::WarningList where warnings are added if the distortion coefficients cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
A std::vector<double> containing the distortion coefficients. The content and size of the vector depend on the distortion model. If the coefficients cannot be parsed, a model-specific default vector is returned, and a warning is added to the list.
-
std::vector<double> getRadialDistortion(nlohmann::json isd, csm::WarningList *list = nullptr)¶
-
std::vector<double> getSunPositions(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the sun positions from an ISD JSON object, converting the values to meters if necessary based on the unit specified in the ISD. The positions are returned as a flat vector of doubles, where each set of three consecutive values represent the x, y, and z coordinates of the sun’s position at a given time.
- Parameters:
isd – The ISD JSON object containing the sun position information and units.
list – Pointer to a csm::WarningList where warnings are added if the sun positions cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
A std::vector<double> containing the sun positions in meters. If the positions cannot be parsed, an empty vector is returned.
-
std::vector<double> getSunVelocities(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the sun velocities from an ISD JSON object, converting the values to meters per second if necessary based on the unit specified in the ISD. The velocities are returned as a flat vector of doubles, where each set of three consecutive values represent the x, y, and z components of the sun’s velocity at a given time.
- Parameters:
isd – The ISD JSON object containing the sun velocity information and units.
list – Pointer to a csm::WarningList where warnings are added if the sun velocities cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
A std::vector<double> containing the sun velocities in meters per second. If the velocities cannot be parsed, an empty vector is returned.
-
std::vector<double> getSensorPositions(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the sensor positions from an ISD JSON object, converting the values to meters if necessary based on the unit specified in the ISD. The positions are returned as a flat vector of doubles, where each set of three consecutive values represent the x, y, and z coordinates of the sensor’s position at a given time.
- Parameters:
isd – The ISD JSON object containing the sensor position information and units.
list – Pointer to a csm::WarningList where warnings are added if the sensor positions cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
A std::vector<double> containing the sensor positions in meters. If the positions cannot be parsed, an empty vector is returned.
-
std::vector<double> getSensorVelocities(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the sensor velocities from an ISD JSON object, converting the values to meters per second if necessary based on the unit specified in the ISD. The velocities are returned as a flat vector of doubles, where each set of three consecutive values represent the x, y, and z components of the sensor’s velocity at a given time.
- Parameters:
isd – The ISD JSON object containing the sensor velocity information and units.
list – Pointer to a csm::WarningList where warnings are added if the sensor velocities cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
A std::vector<double> containing the sensor velocities in meters per second. If the velocities cannot be parsed, an empty vector is returned.
-
std::vector<double> getSensorOrientations(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the sensor orientations from an ISD JSON object. The orientations are returned as a flat vector of doubles, where each set of four consecutive values represent a quaternion (x, y, z, w) describing the sensor’s orientation at a given time.
- Parameters:
isd – The ISD JSON object containing the sensor orientation quaternions.
list – Pointer to a csm::WarningList where warnings are added if the sensor orientations cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
A std::vector<double> containing the sensor orientations as quaternions. If the orientations cannot be parsed, an empty vector is returned.
-
double getWavelength(nlohmann::json isd, csm::WarningList *list = nullptr)¶
@description Retrieves the wavelength from an ISD JSON object. The wavelength is typically specified in meters and represents the central wavelength of the sensor’s operational band.
- Parameters:
isd – The ISD JSON object containing the wavelength.
list – Pointer to a csm::WarningList where warnings are added if the wavelength cannot be parsed from the ISD. This parameter is optional; pass nullptr if you do not need to capture warnings.
- Returns:
The wavelength as a double. If the wavelength cannot be parsed, 0.0 is returned.
-
nlohmann::json stateAsJson(std::string modelState)¶
@description Converts a model state string into a JSON object. This function first sanitizes the input string by removing non-printable characters and then parses the string from the first occurrence of “{” to the last occurrence of “}” into a JSON object.
- Parameters:
modelState – The string representation of the model state to be converted into JSON format.
- Returns:
A JSON object parsed from the sanitized input string. If the string does not contain valid JSON format, the function may throw a json::parse_error.
-
bool readFileInString(std::string const &filename, std::string &str)¶
@description Reads the entire content of a file into a single string. This function clears the output string before reading the file content.
- Parameters:
filename – The path to the file to be read.
str – The output string where the file content will be stored.
- Returns:
A boolean value indicating whether the file was successfully opened and its content read. Returns true if the operation was successful, false otherwise. If the file cannot be opened, an error message is printed to std::cout.
-
void sanitize(std::string &input)¶
@description Sanitizes a given string by replacing characters that are not printable with newlines. This function modifies the original string.
- Parameters:
input – The string to be sanitized.
-
void applyRotationToQuatVec(ale::Rotation const &r, std::vector<double> &quaternions)¶
@description Applies a rotation to a vector of quaternions. The function modifies the quaternions in-place, applying the specified rotation to each quaternion in the vector. Quaternions are expected to be in the format x, y, z, w, but the rotation operation requires them in the format w, x, y, z.
- Parameters:
r – The rotation to be applied, represented as an ale::Rotation object.
quaternions – A reference to a vector of doubles representing the quaternions to be rotated. The vector is modified in-place.
-
void applyRotationTranslationToXyzVec(ale::Rotation const &r, ale::Vec3d const &t, std::vector<double> &xyz)¶
@description Applies a rotation and translation to a vector of xyz vectors. The function modifies the xyz vectors in-place, first applying the specified rotation and then the translation to each vector in the list.
- Parameters:
r – The rotation to be applied, represented as an ale::Rotation object.
t – The translation to be applied, represented as an ale::Vec3d object.
xyz – A reference to a vector of doubles representing the xyz vectors to be transformed. The vector is modified in-place.
-
std::string ephemTimeToCalendarTime(double ephemTime)¶
@description Converts ephemeris time, given in seconds since January 1, 2000 12:00:00 AM GMT, to a calendar time string formatted as “YYYY-MM-DDTHH:MM:SSZ”. The function accounts for the base time and converts the input ephemeris time to the equivalent calendar time in UTC.
- Parameters:
ephemTime – The ephemeris time in seconds since the epoch (January 1, 2000 12:00:00 AM GMT).
- Returns:
A string representing the calendar time in UTC format “YYYY-MM-DDTHH:MM:SSZ”.
-
std::vector<double> pixelToMeter(double line, double sample, std::vector<double> geoTransform)¶
@description Converts pixel coordinates (line, sample) to metric coordinates (x, y) using a given geotransform vector. The geotransform vector contains coefficients that map pixel space to world space.
- Parameters:
line – The line (y) coordinate in pixel space.
sample – The sample (x) coordinate in pixel space.
geoTransform – A vector of doubles representing the geotransformation coefficients.
- Returns:
A vector of doubles containing the metric coordinates (y, x) corresponding to the input pixel coordinates.
-
std::vector<double> meterToPixel(double meter_x, double meter_y, std::vector<double> geoTransform)¶
@description Converts metric coordinates (x, y) to pixel coordinates (line, sample) using a given geotransform vector. The geotransform vector contains coefficients that map world space to pixel space.
- Parameters:
meter_x – The x coordinate in metric space.
meter_y – The y coordinate in metric space.
geoTransform – A vector of doubles representing the geotransformation coefficients.
- Returns:
A vector of doubles containing the pixel coordinates (line, sample) corresponding to the input metric coordinates.