OTB  9.0.0
Orfeo Toolbox
otbLineDetectorImageFilterBase.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 otbLineDetectorImageFilterBase_hxx
22 #define otbLineDetectorImageFilterBase_hxx
23 
25 
26 #include "itkDataObject.h"
27 
28 #include "itkConstantPadImageFilter.h"
29 #include "itkConstNeighborhoodIterator.h"
30 #include "itkNeighborhoodInnerProduct.h"
31 #include "itkImageRegionIterator.h"
32 #include "itkNeighborhoodAlgorithm.h"
33 #include "itkProgressReporter.h"
34 #include "otbMacro.h"
35 #include "otbMath.h"
36 
37 namespace otb
38 {
39 
43 template <class TInputImage, class TOutputImage, class TOutputImageDirection, class InterpolatorType>
45 {
46  this->SetNumberOfRequiredOutputs(2);
47  this->SetNumberOfRequiredOutputs(1);
48  m_Radius.Fill(1);
49  m_LengthLine = 1;
50  m_WidthLine = 0;
51  m_Threshold = 0;
52  m_NumberOfDirections = 4;
53  m_FaceList.Fill(0);
54  // THOMAS : hence base constructor
55  // this->SetNthOutput(0, OutputImageType::New());
56  // this->SetNthOutput(1, OutputImageType::New());
57 }
59 
60 template <class TInputImage, class TOutputImage, class TOutputImageDirection, class InterpolatorType>
62 {
63  // call the superclass' implementation of this method
64  Superclass::GenerateInputRequestedRegion();
65 
66  // get pointers to the input and output
67  typename Superclass::InputImagePointer inputPtr = const_cast<TInputImage*>(this->GetInput());
68  typename Superclass::OutputImagePointer outputPtr = this->GetOutput();
69 
70  if (!inputPtr || !outputPtr)
71  {
72  return;
73  }
74 
75  // get a copy of the input requested region (should equal the output
76  // requested region)
77  typename TInputImage::RegionType inputRequestedRegion;
78  inputRequestedRegion = inputPtr->GetRequestedRegion();
79 
80  // Define the size of the region by the radius
81  m_Radius[1] = static_cast<unsigned int>(3 * (2 * m_WidthLine + 1) + 2);
82  m_Radius[0] = 2 * m_LengthLine + 1;
83 
84  // Define the size of the facelist by taking into account the rotation of the region
85  m_FaceList[0] = static_cast<unsigned int>(std::sqrt(static_cast<double>((m_Radius[0] * m_Radius[0]) + (m_Radius[1] * m_Radius[1]))) + 1);
86  m_FaceList[1] = m_FaceList[0];
87 
88  otbMsgDevMacro(<< "Radius : " << m_Radius[0] << " " << m_Radius[1]);
89  otbMsgDevMacro(<< "-> FaceList : " << m_FaceList[0] << " " << m_FaceList[1]);
90 
91  // pad the input requested region by the operator radius
92  inputRequestedRegion.PadByRadius(m_FaceList);
93 
94  // crop the input requested region at the input's largest possible region
95  if (inputRequestedRegion.Crop(inputPtr->GetLargestPossibleRegion()))
96  {
97  inputPtr->SetRequestedRegion(inputRequestedRegion);
98  return;
99  }
100  else
101  {
102  // Couldn't crop the region (requested region is outside the largest
103  // possible region). Throw an exception.
104 
105  // store what we tried to request (prior to trying to crop)
106  inputPtr->SetRequestedRegion(inputRequestedRegion);
107 
108  // build an exception
109  itk::InvalidRequestedRegionError e(__FILE__, __LINE__);
110  std::ostringstream msg;
111  msg << static_cast<const char*>(this->GetNameOfClass()) << "::GenerateInputRequestedRegion()";
112  e.SetLocation(msg.str());
113  e.SetDescription("Requested region is (at least partially) outside the largest possible region.");
114  e.SetDataObject(inputPtr);
115  throw e;
116  }
117 }
118 
119 /*
120  * Set up state of filter before multi-threading.
121  * InterpolatorType::SetInputImage is not thread-safe and hence
122  * has to be set up before ThreadedGenerateData
123  */
124 template <class TInputImage, class TOutputImage, class TOutputImageDirection, class InterpolatorType>
126 {
127  typename OutputImageType::Pointer output = this->GetOutput();
128  output->FillBuffer(0);
129  typename OutputImageType::Pointer outputDirection = this->GetOutputDirection();
130  outputDirection->FillBuffer(0);
131 }
132 
133 template <class TInputImage, class TOutputImage, class TOutputImageDirection, class InterpolatorType>
135  const OutputImageRegionType& outputRegionForThread, itk::ThreadIdType threadId)
136 {
137 
138  typename InputImageType::ConstPointer input = this->GetInput();
139 
140  InterpolatorPointer interpolator2 = InterpolatorType::New();
141  interpolator2->SetInputImage(input);
142 
143  itk::ZeroFluxNeumannBoundaryCondition<InputImageType> nbc;
144  itk::ConstNeighborhoodIterator<InputImageType> bit, cit;
145  typename itk::ConstNeighborhoodIterator<InputImageType>::OffsetType off;
146  itk::ImageRegionIterator<OutputImageType> it;
147  itk::ImageRegionIterator<OutputImageType> itdir;
148 
149  // Allocate output
150  typename OutputImageType::Pointer output = this->GetOutput();
151  typename OutputImageType::Pointer outputDir = this->GetOutputDirection();
152 
153  // Find the data-set boundary "faces"
154  typename itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<InputImageType>::FaceListType faceList;
155  typename itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<InputImageType>::FaceListType::iterator fit;
156 
157  itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<InputImageType> bC;
158  faceList = bC(input, outputRegionForThread, m_FaceList);
159 
160  // support progress methods/callbacks
161  itk::ProgressReporter progress(this, threadId, outputRegionForThread.GetNumberOfPixels());
162 
163  typename TInputImage::IndexType bitIndex;
164  typename InterpolatorType::ContinuousIndexType Index;
165 
166  // --------------------------------------------------------------------------
167 
168  // Number of direction
169  const unsigned int NB_DIR = this->GetNumberOfDirections();
170  // Number of zone
171  const int NB_ZONE = 3;
172  // Definition of the 4 directions
173  // double Theta[NB_DIR];
174  // ROMAIN
175  double* Theta = new double[NB_DIR];
176 
177  // La rotation nulle correspond a un contour horizontal -> 0 !!
178  for (unsigned int i = 0; i < NB_DIR; ++i)
179  {
180  Theta[i] = (CONST_PI * (i / double(NB_DIR)));
181  /* if(Theta[i]>CONST_PI)
182  Theta[i] = Theta[i]-CONST_PI;
183  if((i/double(NB_DIR))==0.5)
184  Theta[i]=0.; */
185  }
186 
187  // Number of the zone
188  unsigned int zone;
189 
190  // Intensity of the linear feature
191  double R;
192 
193  // Direction of detection
194  double Direction = 0.;
195 
196  // Pixel location in the input image
197  int X, Y;
198 
199  // Pixel location after rotation in the system axis of the region
200  double xout, yout;
201 
202  // location of the central pixel in the input image
203  int Xc, Yc;
204 
205  // location of the central pixel between zone 1 and 2 and between zone 1 and 3
206  int Yc12, Yc13;
207 
208  //---------------------------------------------------------------------------
209  otbMsgDevMacro(<< "Theta : " << Theta[0] << " " << Theta[1] << " " << Theta[2] << " " << Theta[3]);
210 
211  // Process each of the boundary faces. These are N-d regions which border
212  // the edge of the buffer.
213 
214  // bool interiorFace = true;
215 
216  for (fit = faceList.begin(); fit != faceList.end(); ++fit)
217  {
218  bit = itk::ConstNeighborhoodIterator<InputImageType>(m_Radius, input, *fit);
219  cit = itk::ConstNeighborhoodIterator<InputImageType>(m_FaceList, input, *fit);
220 
221  it = itk::ImageRegionIterator<OutputImageType>(output, *fit);
222  itdir = itk::ImageRegionIterator<OutputImageType>(outputDir, *fit);
223 
224  bit.OverrideBoundaryCondition(&nbc);
225  bit.GoToBegin();
226  cit.OverrideBoundaryCondition(&nbc);
227  cit.GoToBegin();
228 
229  otbMsgDevMacro(<< " ------------------- FaceList --------------------------");
230 
231  while ((!bit.IsAtEnd()) && (!cit.IsAtEnd()))
232  {
233  InterpolatorPointer interpolator = InterpolatorType::New();
234  // Location of the central pixel of the region
235  off.Fill(0);
236  bitIndex = bit.GetIndex(off);
237  Xc = bitIndex[0];
238  Yc = bitIndex[1];
239 
240  // JULIEN : If the processed region is the center face
241  // the input image can be used for the interpolation
242  // if (interiorFace)
243  // {
244  // interpolator->SetInputImage(input);
245  // }
246  // else we must feed the interpolator with a partial image corresponding
247  // to the boundary conditions
248  // else
249  // {
250  typename InputImageType::RegionType tempRegion;
251  typename InputImageType::SizeType tempSize;
252  tempSize[0] = 2 * m_FaceList[0] + 1;
253  tempSize[1] = 2 * m_FaceList[1] + 1;
254  tempRegion.SetSize(tempSize);
255  typename itk::ConstNeighborhoodIterator<InputImageType>::OffsetType tempIndex;
256  tempIndex[0] = off[0] - m_FaceList[0];
257  tempIndex[1] = off[1] - m_FaceList[1];
258  tempRegion.SetIndex(cit.GetIndex(tempIndex));
259  typename InputImageType::Pointer tempImage = InputImageType::New();
260  tempImage->SetRegions(tempRegion);
261  tempImage->Allocate();
262 
263  for (unsigned int p = 0; p <= 2 * m_FaceList[0]; ++p)
264  {
265  for (unsigned int q = 0; q <= 2 * m_FaceList[1]; q++)
266  {
267  typename itk::ConstNeighborhoodIterator<InputImageType>::OffsetType index;
268  index[0] = p - m_FaceList[0];
269  index[1] = q - m_FaceList[1];
270  tempImage->SetPixel(cit.GetIndex(index), cit.GetPixel(index));
271  }
272  }
273  interpolator->SetInputImage(tempImage);
274  // }
275 
276  // Location of the central pixel between zone 1 and zone 2
277  Yc12 = Yc - m_WidthLine - 1;
278 
279  // Location of the central pixel between zone 1 and zone 3
280  Yc13 = Yc + m_WidthLine + 1;
281 
282  // Contains for the 4 directions the the pixels belonging to each zone
283  // std::vector<double> PixelValues[NB_DIR][NB_ZONE];
284  // ROMAIN
285  std::vector<double>** PixelValues = nullptr;
286  PixelValues = new std::vector<double>*[NB_DIR];
287  for (unsigned int i = 0; i < NB_DIR; ++i)
288  {
289  PixelValues[i] = nullptr;
290  PixelValues[i] = new std::vector<double>[NB_ZONE];
291  }
292  // otbMsgDevMacro( << "\tCentre Xc/Yc="<<Xc<<" "<<Yc<<" Yc12/Yc13="<<Yc12<<" "<<Yc13);
293  // Loop on the region
294  for (unsigned int i = 0; i < m_Radius[0]; ++i)
295  for (unsigned int j = 0; j < m_Radius[1]; ++j)
296  {
297 
298  off[0] = i - m_Radius[0] / 2;
299  off[1] = j - m_Radius[1] / 2;
300 
301  bitIndex = bit.GetIndex(off);
302  X = bitIndex[0];
303  Y = bitIndex[1];
304 
305  // We determine in the horizontal direction with which zone the pixel belongs.
306  if (Y < Yc12)
307  zone = 1;
308  else if ((Yc12 < Y) && (Y < Yc13))
309  zone = 0;
310  else if (Y > Yc13)
311  zone = 2;
312  else
313  continue;
314  // otbMsgDevMacro( << "\t\tPoint traite (i, j)=("<<i<<","<<j<<") -> X, Y="<<X<<","<<Y<<" zone="<<zone);
315  // Loop on the directions
316  for (unsigned int dir = 0; dir < NB_DIR; ++dir)
317  {
318  // ROTATION( (X-Xc), (Y-Yc), Theta[dir], xout, yout);
319 
320  xout = (X - Xc) * std::cos(Theta[dir]) - (Y - Yc) * std::sin(Theta[dir]);
321  yout = (X - Xc) * std::sin(Theta[dir]) + (Y - Yc) * std::cos(Theta[dir]);
322 
323  Index[0] = static_cast<CoordRepType>(xout + Xc);
324  Index[1] = static_cast<CoordRepType>(yout + Yc);
325 
326  PixelValues[dir][zone].push_back(static_cast<double>(interpolator->EvaluateAtContinuousIndex(Index)));
327  }
328  } // end of the loop on the pixels of the region
329 
330  R = 0.;
331  Direction = 0.;
332 
333  // Loop on the 4 directions
334 
335  for (unsigned int dir = 0; dir < NB_DIR; ++dir)
336  {
337 
338  double Rtemp = this->ComputeMeasure(&PixelValues[dir][0], &PixelValues[dir][1], &PixelValues[dir][2]);
339 
340  if (Rtemp > R)
341  {
342  R = Rtemp;
343  Direction = Theta[dir];
344  }
345 
346  } // end of the loop on the directions
347 
348  // otbMsgDevMacro( << "\t\tR, Direction : "<<R<<","<<Direction);
349  if (R >= this->GetThreshold())
350  {
351 
352  // Assignment of this value to the output pixel
353  it.Set(static_cast<OutputPixelType>(R));
354 
355  // Assignment of this value to the "outputdir" pixel
356  itdir.Set(static_cast<OutputPixelType>(Direction));
357  }
358  else
359  {
360 
361  it.Set(itk::NumericTraits<OutputPixelType>::Zero);
362 
363  itdir.Set(static_cast<OutputPixelType>(0));
364  }
365  ++bit;
366  ++cit;
367  ++it;
368  ++itdir;
369  // interiorFace = false;
370  progress.CompletedPixel();
371 
372  // ROMAIN
373  for (unsigned int i = 0; i < NB_DIR; ++i)
374  {
375  delete[] PixelValues[i];
376  PixelValues[i] = nullptr;
377  }
378  delete[] PixelValues;
379  PixelValues = nullptr;
380  }
381  }
382  delete[] Theta;
383 }
384 
385 template <class TInputImage, class TOutputImage, class TOutputImageDirection, class InterpolatorType>
387  std::vector<double>* itkNotUsed(m2),
388  std::vector<double>* itkNotUsed(m3))
389 {
390  return 0;
391 }
392 
396 template <class TInputImage, class TOutputImage, class TOutputImageDirection, class InterpolatorType>
398 {
399  Superclass::PrintSelf(os, indent);
400  os << indent << "Length: " << m_LengthLine << std::endl;
401  os << indent << "Width: " << m_WidthLine << std::endl;
402 }
404 
405 } // end namespace otb
406 
407 #endif
otb::LineDetectorImageFilterBase::ThreadedGenerateData
void ThreadedGenerateData(const OutputImageRegionType &outputRegionForThread, itk::ThreadIdType threadId) override
Definition: otbLineDetectorImageFilterBase.hxx:134
otb::CONST_PI
constexpr double CONST_PI
Definition: otbMath.h:49
otb::LineDetectorImageFilterBase::ComputeMeasure
virtual double ComputeMeasure(std::vector< double > *m1, std::vector< double > *m2, std::vector< double > *m3)
Definition: otbLineDetectorImageFilterBase.hxx:386
otb::LineDetectorImageFilterBase::LineDetectorImageFilterBase
LineDetectorImageFilterBase()
Definition: otbLineDetectorImageFilterBase.hxx:44
otbMath.h
otb
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
Definition: otbJoinContainer.h:32
otbMacro.h
otb::LineDetectorImageFilterBase::BeforeThreadedGenerateData
void BeforeThreadedGenerateData() override
Definition: otbLineDetectorImageFilterBase.hxx:125
otbLineDetectorImageFilterBase.h
otb::LineDetectorImageFilterBase::GenerateInputRequestedRegion
void GenerateInputRequestedRegion() override
Definition: otbLineDetectorImageFilterBase.hxx:61
otb::LineDetectorImageFilterBase< TInputImage, TOutputImage, TOutputImage, itk::LinearInterpolateImageFunction< TInputImage > >::InterpolatorPointer
InterpolatorType::Pointer InterpolatorPointer
Definition: otbLineDetectorImageFilterBase.h:82
otb::LineDetectorImageFilterBase::PrintSelf
void PrintSelf(std::ostream &os, itk::Indent indent) const override
Definition: otbLineDetectorImageFilterBase.hxx:397
otbMsgDevMacro
#define otbMsgDevMacro(x)
Definition: otbMacro.h:64
otb::LineDetectorImageFilterBase< TInputImage, TOutputImage, TOutputImage, itk::LinearInterpolateImageFunction< TInputImage > >::CoordRepType
InterpolatorType::CoordRepType CoordRepType
Definition: otbLineDetectorImageFilterBase.h:83
otb::LineDetectorImageFilterBase< TInputImage, TOutputImage, TOutputImage, itk::LinearInterpolateImageFunction< TInputImage > >::OutputPixelType
OutputImageType::PixelType OutputPixelType
Definition: otbLineDetectorImageFilterBase.h:89
otb::ImageToModulusAndDirectionImageFilter< TInputImage, TOutputImage, TOutputImage >::OutputImageRegionType
OutputImageType::RegionType OutputImageRegionType
Definition: otbImageToModulusAndDirectionImageFilter.h:71