USGSCSM Astro Frame Sensor Model Class

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

Public Functions

UsgsAstroFrameSensorModel()
~UsgsAstroFrameSensorModel()
bool isValidModelState(const std::string &stringState, csm::WarningList *warnings)
bool isValidIsd(const std::string &stringIsd, csm::WarningList *warnings)
virtual csm::ImageCoord groundToImage(const csm::EcefCoord &groundPt, double desiredPrecision = 0.001, double *achievedPrecision = NULL, csm::WarningList *warnings = NULL) const

UsgsAstroFrameSensorModel::groundToImage.

Parameters
  • groundPt

  • desiredPrecision

  • achievedPrecision

  • warnings

Returns

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

std::string constructStateFromIsd(const std::string &jsonIsd, csm::WarningList *warnings)
void reset()
virtual csm::ImageCoordCovar groundToImage(const csm::EcefCoordCovar &groundPt, double desiredPrecision = 0.001, double *achievedPrecision = NULL, csm::WarningList *warnings = NULL) const
virtual csm::ImageCoord groundToImage(const csm::EcefCoord &ground_pt, const std::vector<double> &adjustments, double desired_precision = 0.001, double *achieved_precision = NULL, csm::WarningList *warnings = NULL) const

UsgsAstroFrameSensorModel::groundToImage.

Parameters
  • groundPt

  • adjustments

  • desired_precision

  • achieved_precision

  • warnings

Returns

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

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

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

Parameters
  • sample – Sample of the input image.

  • line – Line of the input image.

  • height – ???

Returns

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

virtual csm::EcefCoordCovar imageToGround(const csm::ImageCoordCovar &imagePt, double height, double heightVariance, double desiredPrecision = 0.001, double *achievedPrecision = NULL, csm::WarningList *warnings = NULL) const
virtual csm::EcefLocus imageToProximateImagingLocus(const csm::ImageCoord &imagePt, const csm::EcefCoord &groundPt, double desiredPrecision = 0.001, double *achievedPrecision = NULL, csm::WarningList *warnings = NULL) const
virtual csm::EcefLocus imageToRemoteImagingLocus(const csm::ImageCoord &imagePt, double desiredPrecision = 0.001, double *achievedPrecision = NULL, csm::WarningList *warnings = NULL) const
virtual csm::ImageCoord getImageStart() const
virtual csm::ImageVector getImageSize() const
virtual std::pair<csm::ImageCoord, csm::ImageCoord> getValidImageRange() const
virtual std::pair<double, double> getValidHeightRange() const
virtual csm::EcefVector getIlluminationDirection(const csm::EcefCoord &groundPt) const

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

Parameters

groundPt – The ground point to find the illumination vector to.

Returns

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

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

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

Parameters

imagePt – Image coordinate to find the sensor position for.

Throws

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

Returns

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

virtual csm::EcefCoord getSensorPosition(double time) const
virtual csm::EcefVector getSensorVelocity(const csm::ImageCoord &imagePt) const

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

Parameters

imagePt – Image coordinate to find the sensor position for.

Throws

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

Returns

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

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

UsgsAstroFrameSensorModel::computeSensorPartials.

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

Parameters
  • index

  • imagePt

  • groundPt

  • desiredPrecision

  • achievedPrecision

  • warnings

Returns

The partial derivatives in the line,sample directions.

virtual std::vector<csm::RasterGM::SensorPartials> computeAllSensorPartials(const csm::EcefCoord &groundPt, csm::param::Set pSet = csm::param::VALID, double desiredPrecision = 0.001, double *achievedPrecision = NULL, csm::WarningList *warnings = NULL) const
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
virtual std::vector<double> computeGroundPartials(const csm::EcefCoord &groundPt) const
virtual const csm::CorrelationModel &getCorrelationModel() const
virtual std::vector<double> getUnmodeledCrossCovariance(const csm::ImageCoord &pt1, const csm::ImageCoord &pt2) const
virtual csm::Version getVersion() const
virtual std::string getModelName() const
virtual std::string getPedigree() const
virtual std::string getImageIdentifier() const
virtual void setImageIdentifier(const std::string &imageId, csm::WarningList *warnings = NULL)
virtual std::string getSensorIdentifier() const
virtual std::string getPlatformIdentifier() const
virtual std::string getCollectionIdentifier() const
virtual std::string getTrajectoryIdentifier() const
virtual std::string getSensorType() const
virtual std::string getSensorMode() const
virtual std::string getReferenceDateAndTime() const
virtual std::string getModelState() const
virtual void replaceModelState(const std::string &argState)
virtual csm::Ellipsoid getEllipsoid() const
virtual void setEllipsoid(const csm::Ellipsoid &ellipsoid)
virtual csm::EcefCoord getReferencePoint() const
virtual void setReferencePoint(const csm::EcefCoord &groundPt)
virtual int getNumParameters() const
virtual std::string getParameterName(int index) const
virtual std::string getParameterUnits(int index) const
virtual bool hasShareableParameters() const
virtual bool isParameterShareable(int index) const
virtual csm::SharingCriteria getParameterSharingCriteria(int index) const
virtual double getParameterValue(int index) const
virtual void setParameterValue(int index, double value)
virtual csm::param::Type getParameterType(int index) const
virtual void setParameterType(int index, csm::param::Type pType)
virtual double getParameterCovariance(int index1, int index2) const
virtual void setParameterCovariance(int index1, int index2, double covariance)
virtual int getNumGeometricCorrectionSwitches() const
virtual std::string getGeometricCorrectionName(int index) const
virtual void setGeometricCorrectionSwitch(int index, bool value, csm::param::Type pType)
virtual bool getGeometricCorrectionSwitch(int index) const
virtual std::vector<double> getCrossCovarianceMatrix(const GeometricModel &comparisonModel, csm::param::Set pSet = csm::param::VALID, const GeometricModelList &otherModels = GeometricModelList()) const
virtual std::shared_ptr<spdlog::logger> getLogger()
virtual void setLogger(std::string logName)
double getValue(int index, const std::vector<double> &adjustments) const
void calcRotationMatrix(double m[3][3]) const
void calcRotationMatrix(double m[3][3], const std::vector<double> &adjustments) const
void losEllipsoidIntersect(double height, double xc, double yc, double zc, double xl, double yl, double zl, double &x, double &y, double &z, csm::WarningList *warnings = NULL) const

Public Static Functions

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

Public Static Attributes

static const std::string _SENSOR_MODEL_NAME = "USGS_ASTRO_FRAME_SENSOR_MODEL"

Private Members

std::vector<double> m_currentParameterValue
std::vector<double> m_currentParameterCovariance
std::vector<csm::param::Type> m_parameterType
std::vector<double> m_noAdjustments
DistortionType m_distortionType
std::vector<double> m_opticalDistCoeffs
std::vector<double> m_transX
std::vector<double> m_transY
std::vector<double> m_spacecraftVelocity
std::vector<double> m_sunPosition
std::vector<double> m_ccdCenter
std::vector<double> m_iTransS
std::vector<double> m_iTransL
std::vector<double> m_boresight
std::vector<double> m_lineJitter
std::vector<double> m_sampleJitter
std::vector<double> m_lineTimes
double m_majorAxis
double m_minorAxis
double m_focalLength
double m_minElevation
double m_maxElevation
double m_startingDetectorSample
double m_startingDetectorLine
double m_detectorSampleSumming
double m_detectorLineSumming
std::string m_targetName
std::string m_modelName
std::string m_sensorName
std::string m_platformName
std::string m_imageIdentifier
std::string m_collectionIdentifier
double m_ifov
std::string m_instrumentID
double m_focalLengthEpsilon
double m_originalHalfLines
std::string m_spacecraftName
double m_pixelPitch
double m_ephemerisTime
double m_originalHalfSamples
int m_nLines
int m_nSamples
int m_nParameters
csm::EcefCoord m_referencePointXyz
std::shared_ptr<spdlog::logger> m_logger = spdlog::get("usgscsm_logger")
nlohmann::json _state
csm::NoCorrelationModel _no_corr_model

Private Static Attributes

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

Friends

friend class UsgsAstroFramePlugin