OTB  9.0.0
Orfeo Toolbox
otbBandMathImageFilter.hxx
Go to the documentation of this file.
1 /*
2  * Copyright (C) 1999-2011 Insight Software Consortium
3  * Copyright (C) 2005-2022 Centre National d'Etudes Spatiales (CNES)
4  *
5  * This file is part of Orfeo Toolbox
6  *
7  * https://www.orfeo-toolbox.org/
8  *
9  * Licensed under the Apache License, Version 2.0 (the "License");
10  * you may not use this file except in compliance with the License.
11  * You may obtain a copy of the License at
12  *
13  * http://www.apache.org/licenses/LICENSE-2.0
14  *
15  * Unless required by applicable law or agreed to in writing, software
16  * distributed under the License is distributed on an "AS IS" BASIS,
17  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
18  * See the License for the specific language governing permissions and
19  * limitations under the License.
20  */
21 
22 #ifndef otbBandMathImageFilter_hxx
23 #define otbBandMathImageFilter_hxx
24 #include "otbBandMathImageFilter.h"
25 
26 #include "itkImageRegionIterator.h"
27 #include "itkNumericTraits.h"
28 #include "itkProgressReporter.h"
29 #include "otbMacro.h"
30 
31 
32 #include <iostream>
33 #include <string>
34 
35 namespace otb
36 {
37 
39 template <class TImage>
41 {
42  // This number will be incremented each time an image
43  // is added over the one minimumrequired
44  this->SetNumberOfRequiredInputs(1);
45  this->InPlaceOff();
47 
48  m_UnderflowCount = 0;
49  m_OverflowCount = 0;
50  m_ThreadUnderflow.SetSize(1);
51  m_ThreadOverflow.SetSize(1);
52 }
53 
55 template <class TImage>
57 {
58 }
59 
60 template <class TImage>
61 void BandMathImageFilter<TImage>::PrintSelf(std::ostream& os, itk::Indent indent) const
62 {
63  Superclass::PrintSelf(os, indent);
64 
65  os << indent << "Expression: " << m_Expression << std::endl;
66  os << indent << "Computed values follow:" << std::endl;
67  os << indent << "UnderflowCount: " << m_UnderflowCount << std::endl;
68  os << indent << "OverflowCount: " << m_OverflowCount << std::endl;
69  os << indent << "itk::NumericTraits<PixelType>::NonpositiveMin() : " << itk::NumericTraits<PixelType>::NonpositiveMin() << std::endl;
70  os << indent << "itk::NumericTraits<PixelType>::max() : " << itk::NumericTraits<PixelType>::max() << std::endl;
71 }
72 
73 template <class TImage>
75 {
76  this->SetInput(idx, const_cast<TImage*>(image));
77  DataObjectPointerArraySizeType nbInput = this->GetNumberOfInputs();
78  m_VVarName.resize(nbInput + 4);
79  std::ostringstream varName;
80  varName << "b" << nbInput;
81  m_VVarName[idx] = varName.str();
82 
83  m_VVarName[idx + 1] = "idxX";
84  m_VVarName[idx + 2] = "idxY";
85  m_VVarName[idx + 3] = "idxPhyX";
86  m_VVarName[idx + 4] = "idxPhyY";
87 }
88 
89 template <class TImage>
90 void BandMathImageFilter<TImage>::SetNthInput(DataObjectPointerArraySizeType idx, const ImageType* image, const std::string& varName)
91 {
92  this->SetInput(idx, const_cast<TImage*>(image));
93  m_VVarName.resize(this->GetNumberOfInputs() + 4);
94  m_VVarName[idx] = varName;
95 
96  m_VVarName[idx + 1] = "idxX";
97  m_VVarName[idx + 2] = "idxY";
98  m_VVarName[idx + 3] = "idxPhyX";
99  m_VVarName[idx + 4] = "idxPhyY";
100 }
101 
102 template <class TImage>
104 {
105  m_VVarName[idx] = varName;
106 }
107 
108 template <typename TImage>
110 {
111  return const_cast<TImage*>(this->GetInput(idx));
112 }
113 
114 template <typename TImage>
115 void BandMathImageFilter<TImage>::SetExpression(const std::string& expression)
116 {
117  if (m_Expression != expression)
118  m_Expression = expression;
119  this->Modified();
120 }
121 
122 template <typename TImage>
124 {
125  return m_Expression;
126 }
127 
128 template <typename TImage>
130 {
131  return m_VVarName.at(idx);
132 }
133 
134 template <typename TImage>
136 {
137  typename std::vector<ParserType::Pointer>::iterator itParser;
138  unsigned int nbThreads = this->GetNumberOfThreads();
139  unsigned int nbInputImages = this->GetNumberOfInputs();
140  unsigned int nbAccessIndex = 4; // to give access to image and physical index
141  unsigned int i, j;
142  unsigned int inputSize[2];
143  std::vector<std::string> tmpIdxVarNames;
144 
145  tmpIdxVarNames.resize(nbAccessIndex);
146 
147  tmpIdxVarNames.resize(nbAccessIndex);
148  tmpIdxVarNames[0] = "idxX";
149  tmpIdxVarNames[1] = "idxY";
150  tmpIdxVarNames[2] = "idxPhyX";
151  tmpIdxVarNames[3] = "idxPhyY";
152 
153  // Check if input image dimensions matches
154  inputSize[0] = this->GetNthInput(0)->GetLargestPossibleRegion().GetSize(0);
155  inputSize[1] = this->GetNthInput(0)->GetLargestPossibleRegion().GetSize(1);
156 
157  for (unsigned int p = 1; p < nbInputImages; p++)
158  {
159  if ((inputSize[0] != this->GetNthInput(p)->GetLargestPossibleRegion().GetSize(0)) ||
160  (inputSize[1] != this->GetNthInput(p)->GetLargestPossibleRegion().GetSize(1)))
161  {
162  itkExceptionMacro(<< "Input images must have the same dimensions." << std::endl
163  << "band #1 is [" << inputSize[0] << ";" << inputSize[1] << "]" << std::endl
164  << "band #" << p + 1 << " is [" << this->GetNthInput(p)->GetLargestPossibleRegion().GetSize(0) << ";"
165  << this->GetNthInput(p)->GetLargestPossibleRegion().GetSize(1) << "]");
166  }
167  }
168 
169  // Store images specs
170  m_Spacing = this->GetNthInput(0)->GetSignedSpacing();
171  m_Origin = this->GetNthInput(0)->GetOrigin();
172 
173  // Allocate and initialize the thread temporaries
174  m_ThreadUnderflow.SetSize(nbThreads);
175  m_ThreadUnderflow.Fill(0);
176  m_ThreadOverflow.SetSize(nbThreads);
177  m_ThreadOverflow.Fill(0);
178  m_VParser.resize(nbThreads);
179  m_AImage.resize(nbThreads);
180  m_NbVar = nbInputImages + nbAccessIndex;
181  m_VVarName.resize(m_NbVar);
182 
183  for (itParser = m_VParser.begin(); itParser < m_VParser.end(); itParser++)
184  {
185  *itParser = ParserType::New();
186  }
187 
188  for (i = 0; i < nbThreads; ++i)
189  {
190  m_AImage[i].resize(m_NbVar);
191  m_VParser[i]->SetExpr(m_Expression);
192 
193  for (j = 0; j < nbInputImages; ++j)
194  {
195  m_VParser[i]->DefineVar(m_VVarName[j], &(m_AImage[i][j]));
196  }
197 
198  for (j = nbInputImages; j < nbInputImages + nbAccessIndex; ++j)
199  {
200  m_VVarName[j] = tmpIdxVarNames[j - nbInputImages];
201  m_VParser[i]->DefineVar(m_VVarName[j], &(m_AImage[i][j]));
202  }
203  }
204 }
205 
206 template <typename TImage>
208 {
209  unsigned int nbThreads = this->GetNumberOfThreads();
210  unsigned int i;
211 
212  m_UnderflowCount = 0;
213  m_OverflowCount = 0;
214 
215  // Accumulate counts for each thread
216  for (i = 0; i < nbThreads; ++i)
217  {
218  m_UnderflowCount += m_ThreadUnderflow[i];
219  m_OverflowCount += m_ThreadOverflow[i];
220  }
221 
222  if ((m_UnderflowCount != 0) || (m_OverflowCount != 0))
223  otbWarningMacro(<< std::endl
224  << "The Following Parsed Expression : " << this->GetExpression() << std::endl
225  << "Generated " << m_UnderflowCount << " Underflow(s) "
226  << "And " << m_OverflowCount << " Overflow(s) " << std::endl
227  << "The Parsed Expression, The Inputs And The Output "
228  << "Type May Be Incompatible !");
229 }
230 
231 template <typename TImage>
232 void BandMathImageFilter<TImage>::ThreadedGenerateData(const ImageRegionType& outputRegionForThread, itk::ThreadIdType threadId)
233 {
234  double value;
235  unsigned int j;
236  unsigned int nbInputImages = this->GetNumberOfInputs();
237 
238  typedef itk::ImageRegionConstIterator<TImage> ImageRegionConstIteratorType;
239 
240  assert(nbInputImages);
241  std::vector<ImageRegionConstIteratorType> Vit(nbInputImages);
242 
243  for (j = 0; j < nbInputImages; ++j)
244  {
245  Vit[j] = ImageRegionConstIteratorType(this->GetNthInput(j), outputRegionForThread);
246  }
247 
248  itk::ImageRegionIterator<TImage> ot(this->GetOutput(), outputRegionForThread);
249 
250  // support progress methods/callbacks
251  itk::ProgressReporter progress(this, threadId, outputRegionForThread.GetNumberOfPixels());
252 
253  std::vector<double>& threadImage = m_AImage[threadId];
254  ParserType::Pointer const& threadParser = m_VParser[threadId];
255  long& threadUnderflow = m_ThreadUnderflow[threadId];
256  long& threadOverflow = m_ThreadOverflow[threadId];
257  ImageRegionConstIteratorType& firstImageRegion = Vit.front(); // alias for better perfs
258  while (!firstImageRegion.IsAtEnd())
259  {
260  for (j = 0; j < nbInputImages; ++j)
261  {
262  threadImage[j] = static_cast<double>(Vit[j].Get());
263  }
264 
265  // Image Indexes
266  for (j = 0; j < 2; ++j)
267  {
268  threadImage[nbInputImages + j] = static_cast<double>(firstImageRegion.GetIndex()[j]);
269  }
270  for (j = 0; j < 2; ++j)
271  {
272  threadImage[nbInputImages + 2 + j] =
273  static_cast<double>(m_Origin[j]) + static_cast<double>(firstImageRegion.GetIndex()[j]) * static_cast<double>(m_Spacing[j]);
274  }
275 
276  try
277  {
278 
279  value = threadParser->Eval();
280  }
281  catch (itk::ExceptionObject& err)
282  {
283  itkExceptionMacro(<< err);
284  }
285 
286  // Case value is equal to -inf or inferior to the minimum value
287  // allowed by the pixelType cast
288  if (value < double(itk::NumericTraits<PixelType>::NonpositiveMin()))
289  {
290  ot.Set(itk::NumericTraits<PixelType>::NonpositiveMin());
291  threadUnderflow++;
292  }
293  // Case value is equal to inf or superior to the maximum value
294  // allowed by the pixelType cast
295  else if (value > double(itk::NumericTraits<PixelType>::max()))
296  {
297  ot.Set(itk::NumericTraits<PixelType>::max());
298  threadOverflow++;
299  }
300  else
301  {
302  ot.Set(static_cast<PixelType>(value));
303  }
304 
305  for (j = 0; j < nbInputImages; ++j)
306  {
307  ++Vit[j];
308  }
309 
310  ++ot;
311 
312  progress.CompletedPixel();
313  }
314 }
315 
316 } // end namespace otb
317 
318 #endif
otb::BandMathImageFilter::ThreadedGenerateData
void ThreadedGenerateData(const ImageRegionType &outputRegionForThread, itk::ThreadIdType threadId) override
Definition: otbBandMathImageFilter.hxx:232
otb::BandMathImageFilter::GetNthInput
ImageType * GetNthInput(DataObjectPointerArraySizeType idx)
Definition: otbBandMathImageFilter.hxx:109
otb::BandMathImageFilter::SetNthInput
void SetNthInput(DataObjectPointerArraySizeType idx, const ImageType *image)
Definition: otbBandMathImageFilter.hxx:74
otb::BandMathImageFilter::GetExpression
std::string GetExpression() const
Definition: otbBandMathImageFilter.hxx:123
otb::BandMathImageFilter::PrintSelf
void PrintSelf(std::ostream &os, itk::Indent indent) const override
Definition: otbBandMathImageFilter.hxx:61
otbBandMathImageFilter.h
otb::BandMathImageFilter::AfterThreadedGenerateData
void AfterThreadedGenerateData() override
Definition: otbBandMathImageFilter.hxx:207
otb
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Definition: otbJoinContainer.h:32
otb::BandMathImageFilter::ImageRegionType
ImageType::RegionType ImageRegionType
Definition: otbBandMathImageFilter.h:101
otbMacro.h
otb::BandMathImageFilter::GetNthInputName
std::string GetNthInputName(DataObjectPointerArraySizeType idx) const
Definition: otbBandMathImageFilter.hxx:129
otb::BandMathImageFilter::BeforeThreadedGenerateData
void BeforeThreadedGenerateData() override
Definition: otbBandMathImageFilter.hxx:135
otb::BandMathImageFilter::SetExpression
void SetExpression(const std::string &expression)
Definition: otbBandMathImageFilter.hxx:115
otb::BandMathImageFilter::PixelType
ImageType::PixelType PixelType
Definition: otbBandMathImageFilter.h:102
otbWarningMacro
#define otbWarningMacro(x)
Definition: otbMacro.h:65
otb::BandMathImageFilter::DataObjectPointerArraySizeType
itk::ProcessObject::DataObjectPointerArraySizeType DataObjectPointerArraySizeType
Definition: otbBandMathImageFilter.h:107
otb::BandMathImageFilter::ImageType
TImage ImageType
Definition: otbBandMathImageFilter.h:96
otb::BandMathImageFilter::BandMathImageFilter
BandMathImageFilter()
Definition: otbBandMathImageFilter.hxx:40
otb::BandMathImageFilter::SetNthInputName
void SetNthInputName(DataObjectPointerArraySizeType idx, const std::string &expression)
Definition: otbBandMathImageFilter.hxx:103
otb::BandMathImageFilter::~BandMathImageFilter
~BandMathImageFilter() override
Definition: otbBandMathImageFilter.hxx:56
otb::Parser::Pointer
itk::SmartPointer< Self > Pointer
Definition: otbParser.h:50