OTB  10.0.0
Orfeo Toolbox
otbScalarImageToHigherOrderTexturesFilter.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 otbScalarImageToHigherOrderTexturesFilter_hxx
22 #define otbScalarImageToHigherOrderTexturesFilter_hxx
23 
24 #include "otbMacro.h" //for
26 #include "itkImageRegionIteratorWithIndex.h"
27 #include "itkImageRegionIterator.h"
28 #include "itkProgressReporter.h"
29 
30 namespace otb
31 {
32 template <class TInputImage, class TOutputImage>
34  : m_Radius(), m_NumberOfBinsPerAxis(8), m_InputImageMinimum(0), m_InputImageMaximum(255), m_FastCalculations(false), m_SubsampleFactor(), m_SubsampleOffset()
35 {
36  this->DynamicMultiThreadingOn();
37  // There are 10 outputs corresponding to the 8 textures indices
38  this->SetNumberOfRequiredOutputs(10);
39 
40  // Create the 10 outputs
41  this->SetNthOutput(0, OutputImageType::New());
42  this->SetNthOutput(1, OutputImageType::New());
43  this->SetNthOutput(2, OutputImageType::New());
44  this->SetNthOutput(3, OutputImageType::New());
45  this->SetNthOutput(4, OutputImageType::New());
46  this->SetNthOutput(5, OutputImageType::New());
47  this->SetNthOutput(6, OutputImageType::New());
48  this->SetNthOutput(7, OutputImageType::New());
49  this->SetNthOutput(8, OutputImageType::New());
50  this->SetNthOutput(9, OutputImageType::New());
51 
52  m_Radius.Fill(10);
53 
54  // Set the offset directions to their defaults: half of all the possible
55  // directions 1 pixel away. (The other half is included by symmetry.)
56  // We use a neighborhood iterator to calculate the appropriate offsets.
57  typedef itk::Neighborhood<typename InputImageType::PixelType, InputImageType::ImageDimension> NeighborhoodType;
58  NeighborhoodType hood;
59  hood.SetRadius(1);
60 
61  // select all "previous" neighbors that are face+edge+vertex
62  // connected to the current pixel. do not include the center pixel.
63  unsigned int centerIndex = hood.GetCenterNeighborhoodIndex();
64  OffsetVectorPointer offsets = OffsetVector::New();
65  for (unsigned int d = 0; d < centerIndex; d++)
66  {
67  OffsetType offset = hood.GetOffset(d);
68  offsets->push_back(offset);
69  }
70  this->SetOffsets(offsets);
71 
72  this->m_SubsampleFactor.Fill(1);
73  this->m_SubsampleOffset.Fill(0);
74 }
75 
76 template <class TInputImage, class TOutputImage>
78 {
79 }
80 
81 template <class TInputImage, class TOutputImage>
84 {
85  return this->GetOutput(0);
86 }
87 
88 template <class TInputImage, class TOutputImage>
91 {
92  return this->GetOutput(1);
93 }
94 
95 template <class TInputImage, class TOutputImage>
98 {
99  return this->GetOutput(2);
100 }
101 
102 template <class TInputImage, class TOutputImage>
105 {
106  return this->GetOutput(3);
107 }
108 
109 template <class TInputImage, class TOutputImage>
112 {
113  return this->GetOutput(4);
114 }
115 
116 template <class TInputImage, class TOutputImage>
119 {
120  return this->GetOutput(5);
121 }
122 
123 template <class TInputImage, class TOutputImage>
126 {
127  return this->GetOutput(6);
128 }
129 
130 template <class TInputImage, class TOutputImage>
133 {
134  return this->GetOutput(7);
135 }
136 
137 template <class TInputImage, class TOutputImage>
140 {
141  return this->GetOutput(8);
142 }
143 
144 template <class TInputImage, class TOutputImage>
147 {
148  return this->GetOutput(9);
149 }
150 
151 template <class TInputImage, class TOutputImage>
153 {
154  OffsetVectorPointer offsetVector = OffsetVector::New();
155  offsetVector->push_back(offset);
156  this->SetOffsets(offsetVector);
157 }
158 
159 template <class TInputImage, class TOutputImage>
161 {
162  // Compute output size, origin & spacing
163  const InputImageType* inputPtr = this->GetInput();
164  InputRegionType inputRegion = inputPtr->GetLargestPossibleRegion();
165  OutputRegionType outputRegion;
166  outputRegion.SetIndex(0, 0);
167  outputRegion.SetIndex(1, 0);
168  outputRegion.SetSize(0, 1 + (inputRegion.GetSize(0) - 1 - m_SubsampleOffset[0]) / m_SubsampleFactor[0]);
169  outputRegion.SetSize(1, 1 + (inputRegion.GetSize(1) - 1 - m_SubsampleOffset[1]) / m_SubsampleFactor[1]);
170 
171  typename OutputImageType::SpacingType outSpacing = inputPtr->GetSignedSpacing();
172  outSpacing[0] *= m_SubsampleFactor[0];
173  outSpacing[1] *= m_SubsampleFactor[1];
174 
175  typename OutputImageType::PointType outOrigin;
176  inputPtr->TransformIndexToPhysicalPoint(inputRegion.GetIndex() + m_SubsampleOffset, outOrigin);
177 
178  for (unsigned int i = 0; i < this->GetNumberOfOutputs(); i++)
179  {
180  OutputImagePointerType outputPtr = this->GetOutput(i);
181  outputPtr->CopyInformation(inputPtr);
182  outputPtr->SetLargestPossibleRegion(outputRegion);
183  outputPtr->SetOrigin(outOrigin);
184  outputPtr->SetSignedSpacing(outSpacing);
185  }
186 }
187 
188 
189 template <class TInputImage, class TOutputImage>
191 {
192  // First, call superclass implementation
193  Superclass::GenerateInputRequestedRegion();
194 
195  // Retrieve the input and output pointers
196  InputImagePointerType inputPtr = const_cast<InputImageType*>(this->GetInput());
197  OutputImagePointerType outputPtr = this->GetOutput();
198 
199  if (!inputPtr || !outputPtr)
200  {
201  return;
202  }
203 
204  // Retrieve the output requested region
205  // We use only the first output since requested regions for all outputs are enforced to be equal
206  // by the default GenerateOutputRequestedRegiont() implementation
207  OutputRegionType outputRequestedRegion = outputPtr->GetRequestedRegion();
208  typename OutputRegionType::IndexType outputIndex = outputRequestedRegion.GetIndex();
209  typename OutputRegionType::SizeType outputSize = outputRequestedRegion.GetSize();
210  InputRegionType inputLargest = inputPtr->GetLargestPossibleRegion();
211 
212  // Convert index and size to full grid
213  outputIndex[0] = outputIndex[0] * m_SubsampleFactor[0] + m_SubsampleOffset[0] + inputLargest.GetIndex(0);
214  outputIndex[1] = outputIndex[1] * m_SubsampleFactor[1] + m_SubsampleOffset[1] + inputLargest.GetIndex(1);
215  outputSize[0] = 1 + (outputSize[0] - 1) * m_SubsampleFactor[0];
216  outputSize[1] = 1 + (outputSize[1] - 1) * m_SubsampleFactor[1];
217 
218  InputRegionType inputRequestedRegion(outputIndex, outputSize);
219 
220  // Apply the radius
221  inputRequestedRegion.PadByRadius(m_Radius);
222 
223  // Try to apply the requested region to the input image
224  if (inputRequestedRegion.Crop(inputPtr->GetLargestPossibleRegion()))
225  {
226  inputPtr->SetRequestedRegion(inputRequestedRegion);
227  }
228  else
229  {
230  // Build an exception
231  itk::InvalidRequestedRegionError e(__FILE__, __LINE__);
232  e.SetLocation(ITK_LOCATION);
233  e.SetDescription("Requested region is (at least partially) outside the largest possible region.");
234  e.SetDataObject(inputPtr);
235  throw e;
236  }
237 }
238 
239 template <class TInputImage, class TOutputImage>
241 {
242  // Retrieve the input and output pointers
243  const InputImageType* inputPtr = this->GetInput();
244 
245  typedef typename itk::ImageRegionIterator<OutputImageType> IteratorType;
246  std::vector<IteratorType> outputImagesIterators;
247 
248  for (unsigned int i = 0; i < this->GetNumberOfOutputs(); ++i)
249  {
250  outputImagesIterators.push_back(IteratorType(this->GetOutput(i), outputRegionForThread));
251  outputImagesIterators[i].GoToBegin();
252  }
253 
254  // Compute the max possible run length (in physical unit)
255  typename InputImageType::PointType topLeftPoint;
256  typename InputImageType::PointType bottomRightPoint;
257  inputPtr->TransformIndexToPhysicalPoint(outputImagesIterators[0].GetIndex() - m_Radius, topLeftPoint);
258  inputPtr->TransformIndexToPhysicalPoint(outputImagesIterators[0].GetIndex() + m_Radius, bottomRightPoint);
259  double maxDistance = topLeftPoint.EuclideanDistanceTo(bottomRightPoint);
260 
261  InputRegionType inputLargest = inputPtr->GetLargestPossibleRegion();
262 
263  // Iterate on outputs to compute textures
264  while (!outputImagesIterators[0].IsAtEnd())
265  {
266  // Compute the region on which run-length matrix will be estimated
267  typename InputRegionType::IndexType inputIndex;
268  typename InputRegionType::SizeType inputSize;
269 
270  // Convert index to full grid
271  typename OutputImageType::IndexType outIndex;
272 
273  for (unsigned int dim = 0; dim < InputImageType::ImageDimension; ++dim)
274  {
275  outIndex[dim] = outputImagesIterators[0].GetIndex()[dim] * m_SubsampleFactor[dim] + m_SubsampleOffset[dim] + inputLargest.GetIndex(dim);
276  inputIndex[dim] = outIndex[dim] - m_Radius[dim];
277  inputSize[dim] = 2 * m_Radius[dim] + 1;
278  }
279 
280  // Build the input region
281  InputRegionType inputRegion;
282  inputRegion.SetIndex(inputIndex);
283  inputRegion.SetSize(inputSize);
284 
285  inputRegion.Crop(inputPtr->GetBufferedRegion());
286 
287  // Create a local image corresponding to the input region
288  InputImagePointerType localInputImage = InputImageType::New();
289  localInputImage->SetRegions(inputRegion);
290  localInputImage->Allocate();
291  typedef itk::ImageRegionIteratorWithIndex<InputImageType> ImageRegionIteratorType;
292  typedef itk::ImageRegionConstIteratorWithIndex<InputImageType> ImageRegionConstIteratorType;
293  ImageRegionConstIteratorType itInputPtr(inputPtr, inputRegion);
294  ImageRegionIteratorType itLocalInputImage(localInputImage, inputRegion);
295  for (itInputPtr.GoToBegin(), itLocalInputImage.GoToBegin(); !itInputPtr.IsAtEnd(); ++itInputPtr, ++itLocalInputImage)
296  {
297  itLocalInputImage.Set(itInputPtr.Get());
298  }
299 
300  typename ScalarImageToRunLengthFeaturesFilterType::Pointer runLengthFeatureCalculator = ScalarImageToRunLengthFeaturesFilterType::New();
301  runLengthFeatureCalculator->SetInput(localInputImage);
302  runLengthFeatureCalculator->SetOffsets(m_Offsets);
303  runLengthFeatureCalculator->SetNumberOfBinsPerAxis(m_NumberOfBinsPerAxis);
304  runLengthFeatureCalculator->SetPixelValueMinMax(m_InputImageMinimum, m_InputImageMaximum);
305  runLengthFeatureCalculator->SetDistanceValueMinMax(0, maxDistance);
306 
307  runLengthFeatureCalculator->Update();
308 
309  typename ScalarImageToRunLengthFeaturesFilterType::FeatureValueVector& featuresMeans = *(runLengthFeatureCalculator->GetFeatureMeans().GetPointer());
310 
311  // Fill output
312  for (unsigned int i = 0; i < this->GetNumberOfOutputs(); ++i)
313  {
314  // Fill output
315  outputImagesIterators[i].Set(featuresMeans[i]);
316  // Increment iterators
317  ++outputImagesIterators[i];
318  }
319  }
320 }
321 
322 } // End namespace otb
323 
324 #endif
virtual void SetOffsets(const OffsetVector *_arg)
void DynamicThreadedGenerateData(const OutputRegionType &outputRegion) override
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.