OTB  10.0.0
Orfeo Toolbox
Public Types | Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
otb::Spot5SensorModel Class Reference

#include <otbSpot5SensorModel.h>

+ Collaboration diagram for otb::Spot5SensorModel:

Public Types

using ContinuousIndexType = PolygonType::ContinuousIndexType
 
using DurationType = MetaData::Duration
 
using MatrixType = itk::Matrix< double, 3, 3 >
 
using Point2DType = itk::Point< double, 2 >
 
using Point3DType = itk::Point< double, 3 >
 
using PolygonType = otb::Polygon<>
 
using TimeType = MetaData::TimePoint
 
using Vector3DType = itk::Vector< double, 3 >
 

Public Member Functions

MatrixType ComputeSatToOrbRotation (double t) const
 
Point3DType GetAttitude (double time) const
 
Point3DType GetBilinearInterpolation (const double &time, const std::vector< Point3DType > &V, const std::vector< double > &T) const
 
Point3DType GetLagrangeInterpolation (const double &time, const std::vector< Point3DType > &V, const std::vector< double > &T) const
 
void GetPixelLookAngleXY (unsigned int line, double &psiX, double &psiY) const
 
Point3DType GetPositionEcf (double time) const
 
Point3DType GetVelocityEcf (double time) const
 
Ephemeris ImagingRay (Point2DType imPt) const
 
Point3DType IntersectRay (const Ephemeris &imRay) const
 
Point3DType LineSampleHeightToWorld (Point2DType imPt, double heightAboveEllipsoid) const
 
Point3DType LineSampleToWorld (Point2DType imPt) const
 
Point3DType NearestIntersection (const Ephemeris &imRay, double offset) const
 
Spot5SensorModeloperator= (const Spot5SensorModel &)=delete
 
 Spot5SensorModel (const ImageMetadata &imd)
 
 Spot5SensorModel (const Spot5SensorModel &)=delete
 
 Spot5SensorModel (const std::string &productType, const Spot5Param &Spot5Param)
 
void UpdateImageMetadata (ImageMetadata &imd)
 
Point2DType WorldToLineSample (const Point3DType &inGeoPoint) const
 
virtual ~Spot5SensorModel ()=default
 

Private Member Functions

itk::Point< double, 3 > EcefToWorld (const itk::Point< double, 3 > &ecefPoint) const
 
void InitBilinearTransform ()
 
virtual bool insideImage (const itk::Point< double, 2 > point, double epsilon) const
 
ContinuousIndexType Point2DToIndex (const itk::Point< double, 2 > point) const
 
ContinuousIndexType Point3DToIndex (const itk::Point< double, 3 > point) const
 
itk::Point< double, 3 > WorldToEcef (const itk::Point< double, 3 > &worldPoint) const
 

Private Attributes

BilinearProjection::Pointer m_BilinearProj
 
otb::GeocentricTransform< otb::TransformDirection::INVERSE, double >::Pointer m_EcefToWorldTransform
 
PolygonType::Pointer m_GroundRect
 
PolygonType::Pointer m_ImageRect
 
std::string m_ProductType
 
Spot5Param m_Spot5Param
 
otb::GeocentricTransform< otb::TransformDirection::FORWARD, double >::Pointer m_WorldToEcefTransform
 

Static Private Attributes

static constexpr double C = 299792458
 

Detailed Description

Definition at line 39 of file otbSpot5SensorModel.h.

Member Typedef Documentation

◆ ContinuousIndexType

Definition at line 54 of file otbSpot5SensorModel.h.

◆ DurationType

Definition at line 56 of file otbSpot5SensorModel.h.

◆ MatrixType

using otb::Spot5SensorModel::MatrixType = itk::Matrix<double, 3, 3>

Definition at line 51 of file otbSpot5SensorModel.h.

◆ Point2DType

using otb::Spot5SensorModel::Point2DType = itk::Point<double, 2>

SPOT5 sensor class, based on OSSIM Spot5SensorModel with ITK/GDAL implementation.

Definition at line 49 of file otbSpot5SensorModel.h.

◆ Point3DType

using otb::Spot5SensorModel::Point3DType = itk::Point<double, 3>

Definition at line 50 of file otbSpot5SensorModel.h.

◆ PolygonType

Definition at line 53 of file otbSpot5SensorModel.h.

◆ TimeType

Definition at line 55 of file otbSpot5SensorModel.h.

◆ Vector3DType

using otb::Spot5SensorModel::Vector3DType = itk::Vector<double, 3>

Definition at line 52 of file otbSpot5SensorModel.h.

Constructor & Destructor Documentation

◆ Spot5SensorModel() [1/3]

otb::Spot5SensorModel::Spot5SensorModel ( const std::string &  productType,
const Spot5Param Spot5Param 
)

◆ Spot5SensorModel() [2/3]

otb::Spot5SensorModel::Spot5SensorModel ( const ImageMetadata imd)

◆ ~Spot5SensorModel()

virtual otb::Spot5SensorModel::~Spot5SensorModel ( )
virtualdefault

◆ Spot5SensorModel() [3/3]

otb::Spot5SensorModel::Spot5SensorModel ( const Spot5SensorModel )
delete

Member Function Documentation

◆ ComputeSatToOrbRotation()

MatrixType otb::Spot5SensorModel::ComputeSatToOrbRotation ( double  t) const

Compute SatToOrb matrix with an input time.

Parameters
[in]tinput time
Returns
3X3 matrix computed

◆ EcefToWorld()

itk::Point<double, 3> otb::Spot5SensorModel::EcefToWorld ( const itk::Point< double, 3 > &  ecefPoint) const
private

Coordinate transformation from ECEF to geographic.

Parameters
ecefPoint3D position ephemeris coordinates
Returns
itk::Point<double, 3> 3D position wgs84

◆ GetAttitude()

Point3DType otb::Spot5SensorModel::GetAttitude ( double  time) const

Get the Attitude of the sensor at time (interpolation of attitude samples vector from metadata).

Parameters
[in]timeinput time
Returns
3D attitude

◆ GetBilinearInterpolation()

Point3DType otb::Spot5SensorModel::GetBilinearInterpolation ( const double &  time,
const std::vector< Point3DType > &  V,
const std::vector< double > &  T 
) const

Compute the bilinear interpolation from a 3D point vector with a double time vector.

Parameters
[in]timetime to interpolate
[in]V3D point vector
[in]TIndex time double vector
Returns
3D point interpolated

◆ GetLagrangeInterpolation()

Point3DType otb::Spot5SensorModel::GetLagrangeInterpolation ( const double &  time,
const std::vector< Point3DType > &  V,
const std::vector< double > &  T 
) const

Compute the lagrange interpolation from a 3D point vector with a double time vector.

Parameters
[in]timetime to interpolate
[in]V3D point vector
[in]TIndex time double vector
Returns
3D point interpolated

◆ GetPixelLookAngleXY()

void otb::Spot5SensorModel::GetPixelLookAngleXY ( unsigned int  line,
double &  psiX,
double &  psiY 
) const

Get look angles on X and Y axis of the sensor at line (interpolation of angles samples vector from metadata).

Parameters
[in]lineinput line
[out]psiXX angle
[out]psiYY angle

◆ GetPositionEcf()

Point3DType otb::Spot5SensorModel::GetPositionEcf ( double  time) const

Get the 3D position of the sensor at time (interpolation of position samples vector from metadata).

Parameters
[in]timeinput time
Returns
3D position

◆ GetVelocityEcf()

Point3DType otb::Spot5SensorModel::GetVelocityEcf ( double  time) const

Get the 3D velocity of the sensor at time (interpolation of velocity samples vector from metadata).

Parameters
[in]timeinput time
Returns
3D position

◆ ImagingRay()

Ephemeris otb::Spot5SensorModel::ImagingRay ( Point2DType  imPt) const

Compute a ray from sensor position and image point.

Parameters
[in]imPtimage point (col row)
Returns
ephemeris point (3D space point with position and velocity)

◆ InitBilinearTransform()

void otb::Spot5SensorModel::InitBilinearTransform ( )
private

◆ insideImage()

virtual bool otb::Spot5SensorModel::insideImage ( const itk::Point< double, 2 >  point,
double  epsilon 
) const
inlineprivatevirtual

Check if point is inside image.

Parameters
point2D point
epsilonoffset value
Returns
true is inside
false not inside

Definition at line 243 of file otbSpot5SensorModel.h.

References otb::Spot5Param::ImageSize, and m_Spot5Param.

◆ IntersectRay()

Point3DType otb::Spot5SensorModel::IntersectRay ( const Ephemeris imRay) const

Compute world point intersected by image ray.

Parameters
[in]imRayephemeris point (3D space point with position and velocity)
Returns
world point (lat,lon,hgt)

◆ LineSampleHeightToWorld()

Point3DType otb::Spot5SensorModel::LineSampleHeightToWorld ( Point2DType  imPt,
double  heightAboveEllipsoid 
) const

Transform image point (col, row) with a height to world point (lat,lon,hgt).

Parameters
[in]imPtimage point (col, row)
[in]heightAboveEllipsoidheight
Returns
image point and height converted to world point (lat,lon,hgt)

◆ LineSampleToWorld()

Point3DType otb::Spot5SensorModel::LineSampleToWorld ( Point2DType  imPt) const

Transform image point (col, row) to world point (lat,lon,hgt).

Parameters
[in]imPtimage point (col, row)
Returns
image point converted to world point (lat,lon,hgt)

◆ NearestIntersection()

Point3DType otb::Spot5SensorModel::NearestIntersection ( const Ephemeris imRay,
double  offset 
) const

Compute the nearest intersection (world point) between the ray and the ellipsoid.

Parameters
[in]imRayephemeris point (3D space point with position and velocity)
[in]offsetdouble offset
Returns
world point (lat,lon,hgt)

◆ operator=()

Spot5SensorModel& otb::Spot5SensorModel::operator= ( const Spot5SensorModel )
delete

◆ Point2DToIndex()

ContinuousIndexType otb::Spot5SensorModel::Point2DToIndex ( const itk::Point< double, 2 >  point) const
private

Convert 2Dpoint to 2DContinuousIndex.

Parameters
point2Dpoint
Returns
ContinuousIndexType ContinuousIndex

◆ Point3DToIndex()

ContinuousIndexType otb::Spot5SensorModel::Point3DToIndex ( const itk::Point< double, 3 >  point) const
private

Convert 3Dpoint to 2DContinuousIndex.

Parameters
point3Dpoint
Returns
ContinuousIndexType ContinuousIndex

◆ UpdateImageMetadata()

void otb::Spot5SensorModel::UpdateImageMetadata ( ImageMetadata imd)

Update a ImageMetadata object with the stored Spot5Param and GCPs, possibly modified from the original metadata by the Spot5SensorModel.

Parameters
[in]imdImageMetadata to be updated

◆ WorldToEcef()

itk::Point<double, 3> otb::Spot5SensorModel::WorldToEcef ( const itk::Point< double, 3 > &  worldPoint) const
private

Coordinate transformation from geographic to ECEF.

Parameters
worldPoint3D position wgs84
Returns
itk::Point<double, 3> 3D position ephemeris coordinates

◆ WorldToLineSample()

Point2DType otb::Spot5SensorModel::WorldToLineSample ( const Point3DType inGeoPoint) const

Transform world point (lat,lon,hgt) to image point (col,row).

Parameters
[in]inGeoPointground point in WGS84 (lat, lon, hgt) coordinates
Returns
ground point converted to image point (col,row) coordinates

Member Data Documentation

◆ C

constexpr double otb::Spot5SensorModel::C = 299792458
staticconstexprprivate

Definition at line 255 of file otbSpot5SensorModel.h.

◆ m_BilinearProj

BilinearProjection::Pointer otb::Spot5SensorModel::m_BilinearProj
private

Definition at line 258 of file otbSpot5SensorModel.h.

◆ m_EcefToWorldTransform

otb::GeocentricTransform<otb::TransformDirection::INVERSE, double>::Pointer otb::Spot5SensorModel::m_EcefToWorldTransform
private

Definition at line 264 of file otbSpot5SensorModel.h.

◆ m_GroundRect

PolygonType::Pointer otb::Spot5SensorModel::m_GroundRect
private

Definition at line 262 of file otbSpot5SensorModel.h.

◆ m_ImageRect

PolygonType::Pointer otb::Spot5SensorModel::m_ImageRect
private

Definition at line 261 of file otbSpot5SensorModel.h.

◆ m_ProductType

std::string otb::Spot5SensorModel::m_ProductType
private

Definition at line 251 of file otbSpot5SensorModel.h.

◆ m_Spot5Param

Spot5Param otb::Spot5SensorModel::m_Spot5Param
private

Definition at line 252 of file otbSpot5SensorModel.h.

Referenced by insideImage().

◆ m_WorldToEcefTransform

otb::GeocentricTransform<otb::TransformDirection::FORWARD, double>::Pointer otb::Spot5SensorModel::m_WorldToEcefTransform
private

Definition at line 265 of file otbSpot5SensorModel.h.


The documentation for this class was generated from the following file: