OTB  9.0.0
Orfeo Toolbox
otbDEMToImageGenerator.hxx
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2005-2022 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 otbDEMToImageGenerator_hxx
22 #define otbDEMToImageGenerator_hxx
23 
24 #include "otbDEMToImageGenerator.h"
25 #include "otbMacro.h"
26 #include "itkProgressReporter.h"
27 
28 namespace otb
29 {
30 
31 template <class TDEMImage>
33 {
34  m_OutputSpacing[0] = 0.0001;
35  m_OutputSpacing[1] = -0.0001;
36  m_OutputSize[0] = 1;
37  m_OutputSize[1] = 1;
38  m_OutputOrigin[0] = 0;
39  m_OutputOrigin[1] = 0;
40  m_DefaultUnknownValue = itk::NumericTraits<PixelType>::ZeroValue();
41  m_AboveEllipsoid = false;
42 
43  m_Transform = GenericRSTransformType::New();
44 }
45 
46 // GenerateOutputInformation method
47 template <class TDEMImage>
49 {
50  DEMImageType* output = this->GetOutput(0);
51 
52  IndexType start;
53  start[0] = 0;
54  start[1] = 0;
55 
56  // Specify region parameters
57  OutputImageRegionType largestPossibleRegion;
58  largestPossibleRegion.SetSize(m_OutputSize);
59  largestPossibleRegion.SetIndex(start);
60 
61  output->SetLargestPossibleRegion(largestPossibleRegion);
62  output->SetSignedSpacing(m_OutputSpacing);
63  output->SetOrigin(m_OutputOrigin);
64 
65  // Add the metadata set by the user to the output
66  output->m_Imd.Add(MDGeom::ProjectionProj, std::string(m_Transform->GetInputProjectionRef()));
67  if (m_Transform->GetInputImageMetadata() != nullptr)
68  output->m_Imd.Merge(*m_Transform->GetInputImageMetadata());
69 }
70 
71 // InstantiateTransform method
72 template <class TDEMImage>
74 {
75  m_Transform->InstantiateTransform();
76 }
77 
78 template <class TDEMImage>
80 {
81  InstantiateTransform();
82  DEMImagePointerType DEMImage = this->GetOutput();
83 
84  // allocate the output buffer
85  DEMImage->SetBufferedRegion(DEMImage->GetRequestedRegion());
86  DEMImage->Allocate();
87  DEMImage->FillBuffer(0);
88 }
89 
90 
91 template <class TDEMImage>
92 void DEMToImageGenerator<TDEMImage>::ThreadedGenerateData(const OutputImageRegionType& outputRegionForThread, itk::ThreadIdType threadId)
93 {
94  DEMImagePointerType DEMImage = this->GetOutput();
95 
96  // Create an iterator that will walk the output region
97  ImageIteratorType outIt = ImageIteratorType(DEMImage, outputRegionForThread);
98 
99  // support progress methods/callbacks
100  itk::ProgressReporter progress(this, threadId, outputRegionForThread.GetNumberOfPixels());
101 
102  // Walk the output image, evaluating the height at each pixel
103  IndexType currentindex;
104  PointType phyPoint;
105  double height;
106  PointType geoPoint;
107 
108  for (outIt.GoToBegin(); !outIt.IsAtEnd(); ++outIt)
109  {
110  currentindex = outIt.GetIndex();
111  DEMImage->TransformIndexToPhysicalPoint(currentindex, phyPoint);
112 
113 
114  if (m_Transform.IsNotNull())
115  {
116  geoPoint = m_Transform->TransformPoint(phyPoint);
117  if (m_AboveEllipsoid)
118  {
119  height = DEMHandler::GetInstance().GetHeightAboveEllipsoid(geoPoint); // Altitude
120  // calculation
121  }
122  else
123  {
124  height = DEMHandler::GetInstance().GetHeightAboveMSL(geoPoint); // Altitude
125  // calculation
126  }
127  }
128  else
129  {
130  if (m_AboveEllipsoid)
131  {
132  height = DEMHandler::GetInstance().GetHeightAboveEllipsoid(phyPoint); // Altitude
133  // calculation
134  }
135  else
136  {
137  height = DEMHandler::GetInstance().GetHeightAboveMSL(phyPoint); // Altitude
138  // calculation
139  }
140  }
141  /* otbMsgDevMacro(<< "Index : (" << currentindex[0]<< "," << currentindex[1] << ") -> PhyPoint : ("
142  << phyPoint[0] << "," << phyPoint[1] << ") -> GeoPoint: ("
143  << geoPoint[0] << "," << geoPoint[1] << ") -> height" << height);
144  */
145  // otbMsgDevMacro(<< "height" << height);
146  // DEM sets a default value (-32768) at point where it doesn't have altitude information.
147  if (!vnl_math_isnan(height))
148  {
149  // Fill the image
150  DEMImage->SetPixel(currentindex, static_cast<PixelType>(height));
151  }
152  else
153  {
154  // Back to the MNT default value
155  DEMImage->SetPixel(currentindex, m_DefaultUnknownValue);
156  }
157  progress.CompletedPixel();
158  }
159 }
160 
161 template <class TDEMImage>
162 void DEMToImageGenerator<TDEMImage>::PrintSelf(std::ostream& os, itk::Indent indent) const
163 {
164  Superclass::PrintSelf(os, indent);
165 
166  os << indent << "Output Spacing:" << m_OutputSpacing[0] << "," << m_OutputSpacing[1] << std::endl;
167  os << indent << "Output Origin:" << m_OutputOrigin[0] << "," << m_OutputOrigin[1] << std::endl;
168  os << indent << "Output Size:" << m_OutputSize[0] << "," << m_OutputSize[1] << std::endl;
169 }
170 
171 } // namespace otb
172 
173 #endif
otb::DEMToImageGenerator::GenerateOutputInformation
void GenerateOutputInformation() override
Definition: otbDEMToImageGenerator.hxx:48
otb::DEMHandler::GetInstance
static DEMHandler & GetInstance()
otbDEMToImageGenerator.h
otb::DEMToImageGenerator::OutputImageRegionType
Superclass::OutputImageRegionType OutputImageRegionType
Definition: otbDEMToImageGenerator.h:70
otb::DEMHandler::GetHeightAboveMSL
double GetHeightAboveMSL(double lon, double lat) const
otb::DEMToImageGenerator::DEMImagePointerType
DEMImageType::Pointer DEMImagePointerType
Definition: otbDEMToImageGenerator.h:56
otb::MDGeom::ProjectionProj
@ ProjectionProj
otb::DEMToImageGenerator::BeforeThreadedGenerateData
void BeforeThreadedGenerateData() override
Definition: otbDEMToImageGenerator.hxx:79
otb
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Definition: otbJoinContainer.h:32
otbMacro.h
otb::DEMToImageGenerator::DEMToImageGenerator
DEMToImageGenerator()
Definition: otbDEMToImageGenerator.hxx:32
otb::DEMToImageGenerator::DEMImageType
TDEMImage DEMImageType
Definition: otbDEMToImageGenerator.h:55
otb::DEMToImageGenerator::ImageIteratorType
itk::ImageRegionIteratorWithIndex< DEMImageType > ImageIteratorType
Definition: otbDEMToImageGenerator.h:71
otb::DEMToImageGenerator::InstantiateTransform
void InstantiateTransform()
Definition: otbDEMToImageGenerator.hxx:73
otb::DEMToImageGenerator::IndexType
OutputImageType::IndexType IndexType
Definition: otbDEMToImageGenerator.h:69
otb::DEMToImageGenerator::ThreadedGenerateData
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, itk::ThreadIdType threadId) override
Definition: otbDEMToImageGenerator.hxx:92
otb::DEMToImageGenerator::PrintSelf
void PrintSelf(std::ostream &os, itk::Indent indent) const override
Definition: otbDEMToImageGenerator.hxx:162
otb::DEMToImageGenerator::PointType
OutputImageType::PointType PointType
Definition: otbDEMToImageGenerator.h:68
otb::DEMToImageGenerator::PixelType
DEMImageType::PixelType PixelType
Definition: otbDEMToImageGenerator.h:57
otb::DEMHandler::GetHeightAboveEllipsoid
double GetHeightAboveEllipsoid(double lon, double lat) const