21 #ifndef otbDEMToImageGenerator_h
22 #define otbDEMToImageGenerator_h
27 #include "itkImageSource.h"
30 #include "itkImageRegionIteratorWithIndex.h"
50 template <
class TDEMImage>
67 typedef typename OutputImageType::SizeType
SizeType;
68 typedef typename OutputImageType::PointType
PointType;
69 typedef typename OutputImageType::IndexType
IndexType;
95 itkGetConstReferenceMacro(OutputSize,
SizeType);
110 itkSetMacro(AboveEllipsoid,
bool);
111 itkGetMacro(AboveEllipsoid,
bool);
112 itkBooleanMacro(AboveEllipsoid);
115 void InstantiateTransform();
124 m_Transform->SetOutputProjectionRef(ref);
131 return m_Transform->GetOutputProjectionRef();
136 m_Transform->SetInputProjectionRef(ref);
142 return m_Transform->GetInputProjectionRef();
148 return m_Transform->GetOutputImageMetadata();
153 m_Transform->SetOutputImageMetadata(imd);
159 return m_Transform->GetInputImageMetadata();
163 m_Transform->SetInputImageMetadata(imd);
169 template <
class TImageType>
172 this->SetOutputOrigin(image->GetOrigin());
173 this->SetOutputSpacing(image->GetSignedSpacing());
175 this->SetOutputSize(image->GetLargestPossibleRegion().GetSize());
176 this->SetOutputProjectionRef(image->GetProjectionRef());
177 this->SetOutputImageMetadata(&(image->GetImageMetadata()));
180 InstantiateTransform();
189 void PrintSelf(std::ostream& os, itk::Indent indent)
const override;
190 void BeforeThreadedGenerateData()
override;
191 void DynamicThreadedGenerateData(
const OutputImageRegionType& outputRegionForThread)
override;
192 void GenerateOutputInformation()
override;
209 #ifndef OTB_MANUAL_INSTANTIATION
Single access point for DEM data retrieval.
Class to generate an image from DEM data.
OutputImageType::SizeType SizeType
void SetOutputParametersFromImage(const TImageType *image)
~DEMToImageGenerator() override
Superclass::OutputImageRegionType OutputImageRegionType
itk::SmartPointer< Self > Pointer
PixelType m_DefaultUnknownValue
itk::SmartPointer< const Self > ConstPointer
std::string GetInputProjectionRef() const
DEMImageType OutputImageType
GenericRSTransformPointerType m_Transform
otb::DEMHandler DEMHandlerType
DEMToImageGenerator(const Self &)=delete
GenericRSTransform GenericRSTransformType
OutputImageType::IndexType IndexType
void SetOutputImageMetadata(const ImageMetadata *imd)
void SetInputProjectionRef(const std::string &ref)
DEMImageType::PixelType PixelType
void operator=(const Self &)=delete
const ImageMetadata * GetInputImageMetadata() const
void SetInputImageMetadata(const ImageMetadata *imd)
OutputImageType::SpacingType SpacingType
const ImageMetadata * GetOutputImageMetadata() const
OutputImageType::PointType PointType
void SetOutputProjectionRef(const std::string &ref)
itk::ImageSource< DEMImageType > Superclass
GenericRSTransformType::Pointer GenericRSTransformPointerType
Superclass::Pointer OutputImagePointer
itk::ImageRegionIteratorWithIndex< DEMImageType > ImageIteratorType
DEMImageType::Pointer DEMImagePointerType
std::string GetOutputProjectionRef() const
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.