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