OTB  10.0.0
Orfeo Toolbox
otbSpot5SensorModel.h
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2005-2024 Centre National d'Etudes Spatiales (CNES)
3  *
4  * This file is part of Orfeo Toolbox
5  *
6  * https://www.orfeo-toolbox.org/
7  *
8  * Licensed under the Apache License, Version 2.0 (the "License");
9  * you may not use this file except in compliance with the License.
10  * You may obtain a copy of the License at
11  *
12  * http://www.apache.org/licenses/LICENSE-2.0
13  *
14  * Unless required by applicable law or agreed to in writing, software
15  * distributed under the License is distributed on an "AS IS" BASIS,
16  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
17  * See the License for the specific language governing permissions and
18  * limitations under the License.
19  */
20 
21 #ifndef otbSpot5SensorModel_h
22 #define otbSpot5SensorModel_h
23 
24 #include "otbImageMetadata.h"
25 #include "otbSpot5Metadata.h"
26 #include "otbGeocentricTransform.h"
27 #include "otbBilinearProjection.h"
28 #include "otbPolygon.h"
29 
30 
31 #include "itkPoint.h"
32 #include <itkMatrix.h>
33 #include "itkVector.h"
34 #include <itkVersor.h>
35 
36 namespace otb
37 {
38 
40 {
41 
46 public:
47 
48 
49  using Point2DType = itk::Point<double, 2>;
50  using Point3DType = itk::Point<double, 3>;
51  using MatrixType = itk::Matrix<double, 3, 3>;
52  using Vector3DType = itk::Vector<double, 3>;
57 
58  Spot5SensorModel(const std::string & productType, const Spot5Param & Spot5Param);
59 
61  virtual ~Spot5SensorModel() = default;
62 
63  Spot5SensorModel(const Spot5SensorModel&) = delete; // non construction-copyable
64  Spot5SensorModel& operator=(const Spot5SensorModel&) = delete; // non copyable
65 
66 
73  Point2DType WorldToLineSample(const Point3DType& inGeoPoint) const;
74 
83  double heightAboveEllipsoid) const;
84 
92 
100 
108  Point3DType NearestIntersection(const Ephemeris& imRay, double offset) const;
109 
116  Point3DType IntersectRay(const Ephemeris& imRay) const;
117 
127  const std::vector<Point3DType>& V,
128  const std::vector<double>& T) const;
129 
139  const std::vector<Point3DType>& V,
140  const std::vector<double>& T) const;
141 
148  Point3DType GetPositionEcf(double time) const;
149 
156  Point3DType GetVelocityEcf(double time) const;
157 
165  void GetPixelLookAngleXY(unsigned int line, double& psiX, double& psiY) const;
166 
173  Point3DType GetAttitude(double time) const;
174 
182 
183 
191 
192 
193 protected:
194 
195 private:
196 
198 
207  itk::Point<double, 3> EcefToWorld(const itk::Point<double, 3> & ecefPoint) const;
208 
209 
216  itk::Point<double, 3> WorldToEcef(const itk::Point<double, 3> & worldPoint) const;
217 
224  ContinuousIndexType Point2DToIndex(const itk::Point<double, 2> point) const;
225 
232  ContinuousIndexType Point3DToIndex(const itk::Point<double, 3> point) const;
233 
234 
243  virtual bool insideImage(const itk::Point<double, 2> point, double epsilon) const
244  {
245  return (point[0] >= -epsilon) &&
246  (point[0] <= (m_Spot5Param.ImageSize[0]+epsilon-1)) &&
247  (point[1] >= -epsilon) &&
248  (point[1] <= (m_Spot5Param.ImageSize[1]+epsilon-1));
249  }
250 
251  std::string m_ProductType;
253 
254  // Speed of light
255  static constexpr double C = 299792458;
256 
257  // Bilinear transform
259 
260  // Image and ground polygon
263 
266 };
267 
268 }
269 
270 #endif
itk::SmartPointer< Self > Pointer
itk::SmartPointer< Self > Pointer
Generic class containing image metadata used in OTB.
Represents a duration.
Definition: otbDateTime.h:171
Represents a point in Time.
Definition: otbDateTime.h:98
This class represent a 2D polygon.
Definition: otbPolygon.h:45
itk::SmartPointer< Self > Pointer
Definition: otbPolygon.h:50
Superclass::ContinuousIndexType ContinuousIndexType
Definition: otbPolygon.h:63
PolygonType::Pointer m_ImageRect
Point3DType IntersectRay(const Ephemeris &imRay) const
Compute world point intersected by image ray.
Point3DType NearestIntersection(const Ephemeris &imRay, double offset) const
Compute the nearest intersection (world point) between the ray and the ellipsoid.
BilinearProjection::Pointer m_BilinearProj
otb::GeocentricTransform< otb::TransformDirection::FORWARD, double >::Pointer m_WorldToEcefTransform
Point3DType GetPositionEcf(double time) const
Get the 3D position of the sensor at time (interpolation of position samples vector from metadata).
Point3DType GetAttitude(double time) const
Get the Attitude of the sensor at time (interpolation of attitude samples vector from metadata).
Point3DType LineSampleToWorld(Point2DType imPt) const
Transform image point (col, row) to world point (lat,lon,hgt).
itk::Point< double, 2 > Point2DType
SPOT5 sensor class, based on OSSIM Spot5SensorModel with ITK/GDAL implementation.
itk::Point< double, 3 > WorldToEcef(const itk::Point< double, 3 > &worldPoint) const
Coordinate transformation from geographic to ECEF.
ContinuousIndexType Point2DToIndex(const itk::Point< double, 2 > point) const
Convert 2Dpoint to 2DContinuousIndex.
Spot5SensorModel & operator=(const Spot5SensorModel &)=delete
Point3DType 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.
PolygonType::Pointer m_GroundRect
static constexpr double C
Ephemeris ImagingRay(Point2DType imPt) const
Compute a ray from sensor position and image point.
otb::GeocentricTransform< otb::TransformDirection::INVERSE, double >::Pointer m_EcefToWorldTransform
Point3DType GetVelocityEcf(double time) const
Get the 3D velocity of the sensor at time (interpolation of velocity samples vector from metadata).
Point2DType WorldToLineSample(const Point3DType &inGeoPoint) const
Transform world point (lat,lon,hgt) to image point (col,row).
itk::Matrix< double, 3, 3 > MatrixType
itk::Point< double, 3 > EcefToWorld(const itk::Point< double, 3 > &ecefPoint) const
Coordinate transformation from ECEF to geographic.
Point3DType 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.
virtual bool insideImage(const itk::Point< double, 2 > point, double epsilon) const
Check if point is inside image.
itk::Vector< double, 3 > Vector3DType
itk::Point< double, 3 > Point3DType
Spot5SensorModel(const Spot5SensorModel &)=delete
MatrixType ComputeSatToOrbRotation(double t) const
Compute SatToOrb matrix with an input time.
ContinuousIndexType Point3DToIndex(const itk::Point< double, 3 > point) const
Convert 3Dpoint to 2DContinuousIndex.
Point3DType LineSampleHeightToWorld(Point2DType imPt, double heightAboveEllipsoid) const
Transform image point (col, row) with a height to world point (lat,lon,hgt).
Spot5SensorModel(const ImageMetadata &imd)
void 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 me...
PolygonType::ContinuousIndexType ContinuousIndexType
void UpdateImageMetadata(ImageMetadata &imd)
Update a ImageMetadata object with the stored Spot5Param and GCPs, possibly modified from the origina...
virtual ~Spot5SensorModel()=default
Spot5SensorModel(const std::string &productType, const Spot5Param &Spot5Param)
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
This structure is used to handle Ephemeris information.
Spot5 sensors parameters.
Point2DType ImageSize