OTB  9.0.0
Orfeo Toolbox
otbScalarImageToTexturesFilter.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 
22 #ifndef otbScalarImageToTexturesFilter_hxx
23 #define otbScalarImageToTexturesFilter_hxx
24 
26 #include "itkImageRegionIteratorWithIndex.h"
27 #include "itkConstNeighborhoodIterator.h"
28 #include "itkImageRegionIterator.h"
29 #include "itkProgressReporter.h"
30 #include "itkNumericTraits.h"
31 #include <vector>
32 #include <cmath>
33 
34 namespace otb
35 {
36 template <class TInputImage, class TOutputImage>
38  : m_Radius(),
39  m_Offset(),
40  m_NeighborhoodRadius(),
41  m_NumberOfBinsPerAxis(8),
42  m_InputImageMinimum(0),
43  m_InputImageMaximum(255),
44  m_SubsampleFactor(),
45  m_SubsampleOffset()
46 {
47  // There are 8 outputs corresponding to the 8 textures indices
48  this->SetNumberOfRequiredOutputs(8);
49 
50  // Create the 8 outputs
51  this->SetNthOutput(0, OutputImageType::New());
52  this->SetNthOutput(1, OutputImageType::New());
53  this->SetNthOutput(2, OutputImageType::New());
54  this->SetNthOutput(3, OutputImageType::New());
55  this->SetNthOutput(4, OutputImageType::New());
56  this->SetNthOutput(5, OutputImageType::New());
57  this->SetNthOutput(6, OutputImageType::New());
58  this->SetNthOutput(7, OutputImageType::New());
59 
60  this->m_SubsampleFactor.Fill(1);
61  this->m_SubsampleOffset.Fill(0);
62 }
63 
64 template <class TInputImage, class TOutputImage>
66 {
67 }
68 
69 template <class TInputImage, class TOutputImage>
71 {
72  if (this->GetNumberOfOutputs() < 1)
73  {
74  return nullptr;
75  }
76  return static_cast<OutputImageType*>(this->GetOutput(0));
77 }
78 
79 template <class TInputImage, class TOutputImage>
81 {
82  if (this->GetNumberOfOutputs() < 2)
83  {
84  return nullptr;
85  }
86  return static_cast<OutputImageType*>(this->GetOutput(1));
87 }
88 
89 template <class TInputImage, class TOutputImage>
91 {
92  if (this->GetNumberOfOutputs() < 3)
93  {
94  return nullptr;
95  }
96  return static_cast<OutputImageType*>(this->GetOutput(2));
97 }
98 
99 template <class TInputImage, class TOutputImage>
102 {
103  if (this->GetNumberOfOutputs() < 4)
104  {
105  return nullptr;
106  }
107  return static_cast<OutputImageType*>(this->GetOutput(3));
108 }
109 
110 template <class TInputImage, class TOutputImage>
112 {
113  if (this->GetNumberOfOutputs() < 5)
114  {
115  return nullptr;
116  }
117  return static_cast<OutputImageType*>(this->GetOutput(4));
118 }
119 
120 template <class TInputImage, class TOutputImage>
123 {
124  if (this->GetNumberOfOutputs() < 6)
125  {
126  return nullptr;
127  }
128  return static_cast<OutputImageType*>(this->GetOutput(5));
129 }
130 
131 template <class TInputImage, class TOutputImage>
134 {
135  if (this->GetNumberOfOutputs() < 7)
136  {
137  return nullptr;
138  }
139  return static_cast<OutputImageType*>(this->GetOutput(6));
140 }
141 
142 template <class TInputImage, class TOutputImage>
145 {
146  if (this->GetNumberOfOutputs() < 8)
147  {
148  return nullptr;
149  }
150  return static_cast<OutputImageType*>(this->GetOutput(7));
151 }
152 
153 template <class TInputImage, class TOutputImage>
155 {
156  // First, call superclass implementation
157  Superclass::GenerateOutputInformation();
158 
159  // Compute output size, origin & spacing
160  InputRegionType inputRegion = this->GetInput()->GetLargestPossibleRegion();
161  OutputRegionType outputRegion;
162  outputRegion.SetIndex(0, 0);
163  outputRegion.SetIndex(1, 0);
164  outputRegion.SetSize(0, 1 + (inputRegion.GetSize(0) - 1 - m_SubsampleOffset[0]) / m_SubsampleFactor[0]);
165  outputRegion.SetSize(1, 1 + (inputRegion.GetSize(1) - 1 - m_SubsampleOffset[1]) / m_SubsampleFactor[1]);
166 
167  typename OutputImageType::SpacingType outSpacing = this->GetInput()->GetSignedSpacing();
168  outSpacing[0] *= m_SubsampleFactor[0];
169  outSpacing[1] *= m_SubsampleFactor[1];
170 
171  typename OutputImageType::PointType outOrigin;
172  this->GetInput()->TransformIndexToPhysicalPoint(inputRegion.GetIndex() + m_SubsampleOffset, outOrigin);
173 
174  for (unsigned int i = 0; i < this->GetNumberOfOutputs(); i++)
175  {
176  OutputImagePointerType outputPtr = this->GetOutput(i);
177  outputPtr->SetLargestPossibleRegion(outputRegion);
178  outputPtr->SetOrigin(outOrigin);
179  outputPtr->SetSignedSpacing(outSpacing);
180  }
181 }
182 
183 template <class TInputImage, class TOutputImage>
185 {
186  // First, call superclass implementation
187  Superclass::GenerateInputRequestedRegion();
188 
189  // Retrieve the input and output pointers
190  InputImagePointerType inputPtr = const_cast<InputImageType*>(this->GetInput());
191  OutputImagePointerType outputPtr = this->GetOutput();
192 
193  if (!inputPtr || !outputPtr)
194  {
195  return;
196  }
197 
198  // Retrieve the output requested region
199  // We use only the first output since requested regions for all outputs are enforced to be equal
200  // by the default GenerateOutputRequestedRegiont() implementation
201  OutputRegionType outputRequestedRegion = outputPtr->GetRequestedRegion();
202 
203  typename OutputRegionType::IndexType outputIndex = outputRequestedRegion.GetIndex();
204  typename OutputRegionType::SizeType outputSize = outputRequestedRegion.GetSize();
205  typename InputRegionType::IndexType inputIndex;
206  typename InputRegionType::SizeType inputSize;
207  InputRegionType inputLargest = inputPtr->GetLargestPossibleRegion();
208 
209  // Convert index and size to full grid
210  outputIndex[0] = outputIndex[0] * m_SubsampleFactor[0] + m_SubsampleOffset[0] + inputLargest.GetIndex(0);
211  outputIndex[1] = outputIndex[1] * m_SubsampleFactor[1] + m_SubsampleOffset[1] + inputLargest.GetIndex(1);
212  outputSize[0] = 1 + (outputSize[0] - 1) * m_SubsampleFactor[0];
213  outputSize[1] = 1 + (outputSize[1] - 1) * m_SubsampleFactor[1];
214 
215  // First, apply offset
216  for (unsigned int dim = 0; dim < InputImageType::ImageDimension; ++dim)
217  {
218  inputIndex[dim] = std::min(outputIndex[dim], outputIndex[dim] + m_Offset[dim]);
219  inputSize[dim] = std::max(outputIndex[dim] + outputSize[dim], outputIndex[dim] + outputSize[dim] + m_Offset[dim]) - inputIndex[dim];
220  }
221 
222  // Build the input requested region
223  InputRegionType inputRequestedRegion;
224  inputRequestedRegion.SetIndex(inputIndex);
225  inputRequestedRegion.SetSize(inputSize);
226 
227  // Apply the radius
228  inputRequestedRegion.PadByRadius(m_Radius);
229 
230  // Try to apply the requested region to the input image
231  if (inputRequestedRegion.Crop(inputPtr->GetLargestPossibleRegion()))
232  {
233  inputPtr->SetRequestedRegion(inputRequestedRegion);
234  }
235  else
236  {
237  // Build an exception
238  itk::InvalidRequestedRegionError e(__FILE__, __LINE__);
239  e.SetLocation(ITK_LOCATION);
240  e.SetDescription("Requested region is (at least partially) outside the largest possible region.");
241  e.SetDataObject(inputPtr);
242  throw e;
243  }
244 }
245 
246 template <class TInputImage, class TOutputImage>
248 {
249  unsigned int minRadius = 0;
250  for (unsigned int i = 0; i < m_Offset.GetOffsetDimension(); i++)
251  {
252  unsigned int distance = std::abs(m_Offset[i]);
253  if (distance > minRadius)
254  {
255  minRadius = distance;
256  }
257  }
258  m_NeighborhoodRadius.Fill(minRadius);
259 }
260 
261 template <class TInputImage, class TOutputImage>
262 void ScalarImageToTexturesFilter<TInputImage, TOutputImage>::ThreadedGenerateData(const OutputRegionType& outputRegionForThread, itk::ThreadIdType threadId)
263 {
264  // Retrieve the input and output pointers
265  InputImagePointerType inputPtr = const_cast<InputImageType*>(this->GetInput());
266  OutputImagePointerType energyPtr = this->GetEnergyOutput();
267  OutputImagePointerType entropyPtr = this->GetEntropyOutput();
268  OutputImagePointerType correlationPtr = this->GetCorrelationOutput();
269  OutputImagePointerType invDiffMomentPtr = this->GetInverseDifferenceMomentOutput();
270  OutputImagePointerType inertiaPtr = this->GetInertiaOutput();
271  OutputImagePointerType clusterShadePtr = this->GetClusterShadeOutput();
272  OutputImagePointerType clusterProminencePtr = this->GetClusterProminenceOutput();
273  OutputImagePointerType haralickCorPtr = this->GetHaralickCorrelationOutput();
274 
275  // Build output iterators
276  itk::ImageRegionIteratorWithIndex<OutputImageType> energyIt(energyPtr, outputRegionForThread);
277  itk::ImageRegionIterator<OutputImageType> entropyIt(entropyPtr, outputRegionForThread);
278  itk::ImageRegionIterator<OutputImageType> correlationIt(correlationPtr, outputRegionForThread);
279  itk::ImageRegionIterator<OutputImageType> invDiffMomentIt(invDiffMomentPtr, outputRegionForThread);
280  itk::ImageRegionIterator<OutputImageType> inertiaIt(inertiaPtr, outputRegionForThread);
281  itk::ImageRegionIterator<OutputImageType> clusterShadeIt(clusterShadePtr, outputRegionForThread);
282  itk::ImageRegionIterator<OutputImageType> clusterProminenceIt(clusterProminencePtr, outputRegionForThread);
283  itk::ImageRegionIterator<OutputImageType> haralickCorIt(haralickCorPtr, outputRegionForThread);
284 
285  // Go to begin
286  energyIt.GoToBegin();
287  entropyIt.GoToBegin();
288  correlationIt.GoToBegin();
289  invDiffMomentIt.GoToBegin();
290  inertiaIt.GoToBegin();
291  clusterShadeIt.GoToBegin();
292  clusterProminenceIt.GoToBegin();
293  haralickCorIt.GoToBegin();
294 
295  const double log2 = std::log(2.0);
296 
297  InputRegionType inputLargest = inputPtr->GetLargestPossibleRegion();
298 
299  // Set-up progress reporting
300  itk::ProgressReporter progress(this, threadId, outputRegionForThread.GetNumberOfPixels());
301 
302  // Iterate on outputs to compute textures
303  while (!energyIt.IsAtEnd() && !entropyIt.IsAtEnd() && !correlationIt.IsAtEnd() && !invDiffMomentIt.IsAtEnd() && !inertiaIt.IsAtEnd() &&
304  !clusterShadeIt.IsAtEnd() && !clusterProminenceIt.IsAtEnd() && !haralickCorIt.IsAtEnd())
305  {
306  // Compute the region on which co-occurence will be estimated
307  typename InputRegionType::IndexType inputIndex;
308  typename InputRegionType::SizeType inputSize;
309 
310  // Convert index to full grid
311  typename OutputImageType::IndexType outIndex;
312 
313  // First, create an window for neighborhood iterator based on m_Radius
314  // For example, if xradius and yradius is 2. window size is 5x5 (2 * radius + 1).
315  for (unsigned int dim = 0; dim < InputImageType::ImageDimension; ++dim)
316  {
317  outIndex[dim] = energyIt.GetIndex()[dim] * m_SubsampleFactor[dim] + m_SubsampleOffset[dim] + inputLargest.GetIndex(dim);
318  inputIndex[dim] = outIndex[dim] - m_Radius[dim];
319  inputSize[dim] = 2 * m_Radius[dim] + 1;
320  }
321 
322  // Build the input region
323  InputRegionType inputRegion;
324  inputRegion.SetIndex(inputIndex);
325  inputRegion.SetSize(inputSize);
326  inputRegion.Crop(inputPtr->GetRequestedRegion());
327 
328  CooccurrenceIndexedListPointerType GLCIList = CooccurrenceIndexedListType::New();
329  GLCIList->Initialize(m_NumberOfBinsPerAxis, m_InputImageMinimum, m_InputImageMaximum);
330 
331  typedef itk::ConstNeighborhoodIterator<InputImageType> NeighborhoodIteratorType;
332  NeighborhoodIteratorType neighborIt;
333  neighborIt = NeighborhoodIteratorType(m_NeighborhoodRadius, inputPtr, inputRegion);
334  for (neighborIt.GoToBegin(); !neighborIt.IsAtEnd(); ++neighborIt)
335  {
336  const InputPixelType centerPixelIntensity = neighborIt.GetCenterPixel();
337  bool pixelInBounds;
338  const InputPixelType pixelIntensity = neighborIt.GetPixel(m_Offset, pixelInBounds);
339  if (!pixelInBounds)
340  {
341  continue; // don't put a pixel in the co-occurrence list if the value is
342  // out of bounds
343  }
344  GLCIList->AddPixelPair(centerPixelIntensity, pixelIntensity);
345  }
346 
347  double pixelMean = 0.;
348  double marginalMean;
349  double marginalDevSquared = 0.;
350  double pixelVariance = 0.;
351 
352  // Create and Initialize marginalSums
353  std::vector<double> marginalSums(m_NumberOfBinsPerAxis, 0);
354 
355  // get co-occurrence vector and totalfrequency
356  VectorType glcVector = GLCIList->GetVector();
357  double totalFrequency = static_cast<double>(GLCIList->GetTotalFrequency());
358 
359  // Normalize the co-occurrence indexed list and compute mean, marginalSum
360  typename VectorType::iterator it = glcVector.begin();
361  while (it != glcVector.end())
362  {
363  double frequency = (*it).second / totalFrequency;
364  CooccurrenceIndexType index = (*it).first;
365  pixelMean += index[0] * frequency;
366  marginalSums[index[0]] += frequency;
367  ++it;
368  }
369 
370  /* Now get the mean and deviaton of the marginal sums.
371  Compute incremental mean and SD, a la Knuth, "The Art of Computer
372  Programming, Volume 2: Seminumerical Algorithms", section 4.2.2.
373  Compute mean and standard deviation using the recurrence relation:
374  M(1) = x(1), M(k) = M(k-1) + (x(k) - M(k-1) ) / k
375  S(1) = 0, S(k) = S(k-1) + (x(k) - M(k-1)) * (x(k) - M(k))
376  for 2 <= k <= n, then
377  sigma = std::sqrt(S(n) / n) (or divide by n-1 for sample SD instead of
378  population SD).
379  */
380  std::vector<double>::const_iterator msIt = marginalSums.begin();
381  marginalMean = *msIt;
382  // Increment iterator to start with index 1
383  ++msIt;
384  for (int k = 2; msIt != marginalSums.end(); ++k, ++msIt)
385  {
386  double M_k_minus_1 = marginalMean;
387  double S_k_minus_1 = marginalDevSquared;
388  double x_k = *msIt;
389  double M_k = M_k_minus_1 + (x_k - M_k_minus_1) / k;
390  double S_k = S_k_minus_1 + (x_k - M_k_minus_1) * (x_k - M_k);
391  marginalMean = M_k;
392  marginalDevSquared = S_k;
393  }
394  marginalDevSquared = marginalDevSquared / m_NumberOfBinsPerAxis;
395 
396  VectorConstIteratorType constVectorIt;
397  constVectorIt = glcVector.begin();
398  while (constVectorIt != glcVector.end())
399  {
400  RelativeFrequencyType frequency = (*constVectorIt).second / totalFrequency;
401  CooccurrenceIndexType index = (*constVectorIt).first;
402  pixelVariance += (index[0] - pixelMean) * (index[0] - pixelMean) * frequency;
403  ++constVectorIt;
404  }
405 
406  double pixelVarianceSquared = pixelVariance * pixelVariance;
407  // Variance is only used in correlation. If variance is 0, then (index[0] - pixelMean) * (index[1] - pixelMean)
408  // should be zero as well. In this case, set the variance to 1. in order to
409  // avoid NaN correlation.
410  if (pixelVarianceSquared < GetPixelValueTolerance())
411  {
412  pixelVarianceSquared = 1.;
413  }
414 
415  // Initialize texture variables;
416  PixelValueType energy = itk::NumericTraits<PixelValueType>::Zero;
417  PixelValueType entropy = itk::NumericTraits<PixelValueType>::Zero;
418  PixelValueType correlation = itk::NumericTraits<PixelValueType>::Zero;
419  PixelValueType inverseDifferenceMoment = itk::NumericTraits<PixelValueType>::Zero;
420  PixelValueType inertia = itk::NumericTraits<PixelValueType>::Zero;
421  PixelValueType clusterShade = itk::NumericTraits<PixelValueType>::Zero;
422  PixelValueType clusterProminence = itk::NumericTraits<PixelValueType>::Zero;
423  PixelValueType haralickCorrelation = itk::NumericTraits<PixelValueType>::Zero;
424 
425  // Compute textures
426  constVectorIt = glcVector.begin();
427  while (constVectorIt != glcVector.end())
428  {
429  CooccurrenceIndexType index = (*constVectorIt).first;
430  RelativeFrequencyType frequency = (*constVectorIt).second / totalFrequency;
431  energy += frequency * frequency;
432  entropy -= (frequency > GetPixelValueTolerance()) ? frequency * std::log(frequency) / log2 : 0;
433  correlation += ((index[0] - pixelMean) * (index[1] - pixelMean) * frequency) / pixelVarianceSquared;
434  inverseDifferenceMoment += frequency / (1.0 + (index[0] - index[1]) * (index[0] - index[1]));
435  inertia += (index[0] - index[1]) * (index[0] - index[1]) * frequency;
436  clusterShade += std::pow((index[0] - pixelMean) + (index[1] - pixelMean), 3) * frequency;
437  clusterProminence += std::pow((index[0] - pixelMean) + (index[1] - pixelMean), 4) * frequency;
438  haralickCorrelation += index[0] * index[1] * frequency;
439  ++constVectorIt;
440  }
441 
442  haralickCorrelation = (fabs(marginalDevSquared) > 1E-8) ? (haralickCorrelation - marginalMean * marginalMean) / marginalDevSquared : 0;
443 
444  // Fill outputs
445  energyIt.Set(energy);
446  entropyIt.Set(entropy);
447  correlationIt.Set(correlation);
448  invDiffMomentIt.Set(inverseDifferenceMoment);
449  inertiaIt.Set(inertia);
450  clusterShadeIt.Set(clusterShade);
451  clusterProminenceIt.Set(clusterProminence);
452  haralickCorIt.Set(haralickCorrelation);
453 
454  // Update progress
455  progress.CompletedPixel();
456 
457  // Increment iterators
458  ++energyIt;
459  ++entropyIt;
460  ++correlationIt;
461  ++invDiffMomentIt;
462  ++inertiaIt;
463  ++clusterShadeIt;
464  ++clusterProminenceIt;
465  ++haralickCorIt;
466  }
467 }
468 
469 } // End namespace otb
470 
471 #endif
otb::ScalarImageToTexturesFilter::GetEnergyOutput
OutputImageType * GetEnergyOutput()
Definition: otbScalarImageToTexturesFilter.hxx:70
otb::ScalarImageToTexturesFilter::GetClusterShadeOutput
OutputImageType * GetClusterShadeOutput()
Definition: otbScalarImageToTexturesFilter.hxx:122
otb::ScalarImageToTexturesFilter::m_SubsampleOffset
OffsetType m_SubsampleOffset
Definition: otbScalarImageToTexturesFilter.h:263
otb
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Definition: otbJoinContainer.h:32
otb::ScalarImageToTexturesFilter::OutputImagePointerType
OutputImageType::Pointer OutputImagePointerType
Definition: otbScalarImageToTexturesFilter.h:129
otb::ScalarImageToTexturesFilter::BeforeThreadedGenerateData
void BeforeThreadedGenerateData() override
Definition: otbScalarImageToTexturesFilter.hxx:247
otb::ScalarImageToTexturesFilter::GetEntropyOutput
OutputImageType * GetEntropyOutput()
Definition: otbScalarImageToTexturesFilter.hxx:80
otb::ScalarImageToTexturesFilter::GetInverseDifferenceMomentOutput
OutputImageType * GetInverseDifferenceMomentOutput()
Definition: otbScalarImageToTexturesFilter.hxx:101
otb::ScalarImageToTexturesFilter::VectorConstIteratorType
VectorType::const_iterator VectorConstIteratorType
Definition: otbScalarImageToTexturesFilter.h:141
otb::ScalarImageToTexturesFilter::m_SubsampleFactor
SizeType m_SubsampleFactor
Definition: otbScalarImageToTexturesFilter.h:260
otb::ScalarImageToTexturesFilter::CooccurrenceIndexedListPointerType
CooccurrenceIndexedListType::Pointer CooccurrenceIndexedListPointerType
Definition: otbScalarImageToTexturesFilter.h:133
otb::ScalarImageToTexturesFilter::VectorType
CooccurrenceIndexedListType::VectorType VectorType
Definition: otbScalarImageToTexturesFilter.h:138
otb::ScalarImageToTexturesFilter::~ScalarImageToTexturesFilter
~ScalarImageToTexturesFilter() override
Definition: otbScalarImageToTexturesFilter.hxx:65
otb::ScalarImageToTexturesFilter::GenerateOutputInformation
void GenerateOutputInformation() override
Definition: otbScalarImageToTexturesFilter.hxx:154
otb::ScalarImageToTexturesFilter::InputRegionType
InputImageType::RegionType InputRegionType
Definition: otbScalarImageToTexturesFilter.h:124
otb::ScalarImageToTexturesFilter::ThreadedGenerateData
void ThreadedGenerateData(const OutputRegionType &outputRegion, itk::ThreadIdType threadId) override
Definition: otbScalarImageToTexturesFilter.hxx:262
otb::ScalarImageToTexturesFilter::GetInertiaOutput
OutputImageType * GetInertiaOutput()
Definition: otbScalarImageToTexturesFilter.hxx:111
otb::ScalarImageToTexturesFilter::InputImagePointerType
InputImageType::Pointer InputImagePointerType
Definition: otbScalarImageToTexturesFilter.h:122
otb::ScalarImageToTexturesFilter::OutputRegionType
OutputImageType::RegionType OutputRegionType
Definition: otbScalarImageToTexturesFilter.h:130
otb::ScalarImageToTexturesFilter::RelativeFrequencyType
CooccurrenceIndexedListType::RelativeFrequencyType RelativeFrequencyType
Definition: otbScalarImageToTexturesFilter.h:137
otb::ScalarImageToTexturesFilter::GetClusterProminenceOutput
OutputImageType * GetClusterProminenceOutput()
Definition: otbScalarImageToTexturesFilter.hxx:133
otb::ScalarImageToTexturesFilter::CooccurrenceIndexType
CooccurrenceIndexedListType::IndexType CooccurrenceIndexType
Definition: otbScalarImageToTexturesFilter.h:135
otb::ScalarImageToTexturesFilter::PixelValueType
CooccurrenceIndexedListType::PixelValueType PixelValueType
Definition: otbScalarImageToTexturesFilter.h:136
otbScalarImageToTexturesFilter.h
otb::ScalarImageToTexturesFilter::GetHaralickCorrelationOutput
OutputImageType * GetHaralickCorrelationOutput()
Definition: otbScalarImageToTexturesFilter.hxx:144
otb::ScalarImageToTexturesFilter::OutputImageType
TOutputImage OutputImageType
Definition: otbScalarImageToTexturesFilter.h:128
otb::ScalarImageToTexturesFilter::InputImageType
TInpuImage InputImageType
Definition: otbScalarImageToTexturesFilter.h:118
otb::ScalarImageToTexturesFilter::InputPixelType
InputImageType::PixelType InputPixelType
Definition: otbScalarImageToTexturesFilter.h:123
otb::ScalarImageToTexturesFilter::ScalarImageToTexturesFilter
ScalarImageToTexturesFilter()
Definition: otbScalarImageToTexturesFilter.hxx:37
otb::ScalarImageToTexturesFilter::GetCorrelationOutput
OutputImageType * GetCorrelationOutput()
Definition: otbScalarImageToTexturesFilter.hxx:90
otb::ScalarImageToTexturesFilter::GenerateInputRequestedRegion
void GenerateInputRequestedRegion() override
Definition: otbScalarImageToTexturesFilter.hxx:184