OTB  10.0.0
Orfeo Toolbox
otbLineDetectorImageFilterBase.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 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)
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 
161  typename TInputImage::IndexType bitIndex;
162  typename InterpolatorType::ContinuousIndexType Index;
163 
164  // --------------------------------------------------------------------------
165 
166  // Number of direction
167  const unsigned int NB_DIR = this->GetNumberOfDirections();
168  // Number of zone
169  const int NB_ZONE = 3;
170  // Definition of the 4 directions
171  // double Theta[NB_DIR];
172  // ROMAIN
173  double* Theta = new double[NB_DIR];
174 
175  // La rotation nulle correspond a un contour horizontal -> 0 !!
176  for (unsigned int i = 0; i < NB_DIR; ++i)
177  {
178  Theta[i] = (CONST_PI * (i / double(NB_DIR)));
179  /* if(Theta[i]>CONST_PI)
180  Theta[i] = Theta[i]-CONST_PI;
181  if((i/double(NB_DIR))==0.5)
182  Theta[i]=0.; */
183  }
184 
185  // Number of the zone
186  unsigned int zone;
187 
188  // Intensity of the linear feature
189  double R;
190 
191  // Direction of detection
192  double Direction = 0.;
193 
194  // Pixel location in the input image
195  int X, Y;
196 
197  // Pixel location after rotation in the system axis of the region
198  double xout, yout;
199 
200  // location of the central pixel in the input image
201  int Xc, Yc;
202 
203  // location of the central pixel between zone 1 and 2 and between zone 1 and 3
204  int Yc12, Yc13;
205 
206  //---------------------------------------------------------------------------
207  otbMsgDevMacro(<< "Theta : " << Theta[0] << " " << Theta[1] << " " << Theta[2] << " " << Theta[3]);
208 
209  // Process each of the boundary faces. These are N-d regions which border
210  // the edge of the buffer.
211 
212  // bool interiorFace = true;
213 
214  for (fit = faceList.begin(); fit != faceList.end(); ++fit)
215  {
216  bit = itk::ConstNeighborhoodIterator<InputImageType>(m_Radius, input, *fit);
217  cit = itk::ConstNeighborhoodIterator<InputImageType>(m_FaceList, input, *fit);
218 
219  it = itk::ImageRegionIterator<OutputImageType>(output, *fit);
220  itdir = itk::ImageRegionIterator<OutputImageType>(outputDir, *fit);
221 
222  bit.OverrideBoundaryCondition(&nbc);
223  bit.GoToBegin();
224  cit.OverrideBoundaryCondition(&nbc);
225  cit.GoToBegin();
226 
227  otbMsgDevMacro(<< " ------------------- FaceList --------------------------");
228 
229  while ((!bit.IsAtEnd()) && (!cit.IsAtEnd()))
230  {
231  InterpolatorPointer interpolator = InterpolatorType::New();
232  // Location of the central pixel of the region
233  off.Fill(0);
234  bitIndex = bit.GetIndex(off);
235  Xc = bitIndex[0];
236  Yc = bitIndex[1];
237 
238  // JULIEN : If the processed region is the center face
239  // the input image can be used for the interpolation
240  // if (interiorFace)
241  // {
242  // interpolator->SetInputImage(input);
243  // }
244  // else we must feed the interpolator with a partial image corresponding
245  // to the boundary conditions
246  // else
247  // {
248  typename InputImageType::RegionType tempRegion;
249  typename InputImageType::SizeType tempSize;
250  tempSize[0] = 2 * m_FaceList[0] + 1;
251  tempSize[1] = 2 * m_FaceList[1] + 1;
252  tempRegion.SetSize(tempSize);
253  typename itk::ConstNeighborhoodIterator<InputImageType>::OffsetType tempIndex;
254  tempIndex[0] = off[0] - m_FaceList[0];
255  tempIndex[1] = off[1] - m_FaceList[1];
256  tempRegion.SetIndex(cit.GetIndex(tempIndex));
257  typename InputImageType::Pointer tempImage = InputImageType::New();
258  tempImage->SetRegions(tempRegion);
259  tempImage->Allocate();
260 
261  for (unsigned int p = 0; p <= 2 * m_FaceList[0]; ++p)
262  {
263  for (unsigned int q = 0; q <= 2 * m_FaceList[1]; q++)
264  {
265  typename itk::ConstNeighborhoodIterator<InputImageType>::OffsetType index;
266  index[0] = p - m_FaceList[0];
267  index[1] = q - m_FaceList[1];
268  tempImage->SetPixel(cit.GetIndex(index), cit.GetPixel(index));
269  }
270  }
271  interpolator->SetInputImage(tempImage);
272  // }
273 
274  // Location of the central pixel between zone 1 and zone 2
275  Yc12 = Yc - m_WidthLine - 1;
276 
277  // Location of the central pixel between zone 1 and zone 3
278  Yc13 = Yc + m_WidthLine + 1;
279 
280  // Contains for the 4 directions the the pixels belonging to each zone
281  // std::vector<double> PixelValues[NB_DIR][NB_ZONE];
282  // ROMAIN
283  std::vector<double>** PixelValues = nullptr;
284  PixelValues = new std::vector<double>*[NB_DIR];
285  for (unsigned int i = 0; i < NB_DIR; ++i)
286  {
287  PixelValues[i] = nullptr;
288  PixelValues[i] = new std::vector<double>[NB_ZONE];
289  }
290  // otbMsgDevMacro( << "\tCentre Xc/Yc="<<Xc<<" "<<Yc<<" Yc12/Yc13="<<Yc12<<" "<<Yc13);
291  // Loop on the region
292  for (unsigned int i = 0; i < m_Radius[0]; ++i)
293  for (unsigned int j = 0; j < m_Radius[1]; ++j)
294  {
295 
296  off[0] = i - m_Radius[0] / 2;
297  off[1] = j - m_Radius[1] / 2;
298 
299  bitIndex = bit.GetIndex(off);
300  X = bitIndex[0];
301  Y = bitIndex[1];
302 
303  // We determine in the horizontal direction with which zone the pixel belongs.
304  if (Y < Yc12)
305  zone = 1;
306  else if ((Yc12 < Y) && (Y < Yc13))
307  zone = 0;
308  else if (Y > Yc13)
309  zone = 2;
310  else
311  continue;
312  // otbMsgDevMacro( << "\t\tPoint traite (i, j)=("<<i<<","<<j<<") -> X, Y="<<X<<","<<Y<<" zone="<<zone);
313  // Loop on the directions
314  for (unsigned int dir = 0; dir < NB_DIR; ++dir)
315  {
316  // ROTATION( (X-Xc), (Y-Yc), Theta[dir], xout, yout);
317 
318  xout = (X - Xc) * std::cos(Theta[dir]) - (Y - Yc) * std::sin(Theta[dir]);
319  yout = (X - Xc) * std::sin(Theta[dir]) + (Y - Yc) * std::cos(Theta[dir]);
320 
321  Index[0] = static_cast<CoordRepType>(xout + Xc);
322  Index[1] = static_cast<CoordRepType>(yout + Yc);
323 
324  PixelValues[dir][zone].push_back(static_cast<double>(interpolator->EvaluateAtContinuousIndex(Index)));
325  }
326  } // end of the loop on the pixels of the region
327 
328  R = 0.;
329  Direction = 0.;
330 
331  // Loop on the 4 directions
332 
333  for (unsigned int dir = 0; dir < NB_DIR; ++dir)
334  {
335 
336  double Rtemp = this->ComputeMeasure(&PixelValues[dir][0], &PixelValues[dir][1], &PixelValues[dir][2]);
337 
338  if (Rtemp > R)
339  {
340  R = Rtemp;
341  Direction = Theta[dir];
342  }
343 
344  } // end of the loop on the directions
345 
346  // otbMsgDevMacro( << "\t\tR, Direction : "<<R<<","<<Direction);
347  if (R >= this->GetThreshold())
348  {
349 
350  // Assignment of this value to the output pixel
351  it.Set(static_cast<OutputPixelType>(R));
352 
353  // Assignment of this value to the "outputdir" pixel
354  itdir.Set(static_cast<OutputPixelType>(Direction));
355  }
356  else
357  {
358 
359  it.Set(itk::NumericTraits<OutputPixelType>::Zero);
360 
361  itdir.Set(static_cast<OutputPixelType>(0));
362  }
363  ++bit;
364  ++cit;
365  ++it;
366  ++itdir;
367  // interiorFace = false;
368 
369  for (unsigned int i = 0; i < NB_DIR; ++i)
370  {
371  delete[] PixelValues[i];
372  PixelValues[i] = nullptr;
373  }
374  delete[] PixelValues;
375  PixelValues = nullptr;
376  }
377  }
378  delete[] Theta;
379 }
380 
381 template <class TInputImage, class TOutputImage, class TOutputImageDirection, class InterpolatorType>
383  std::vector<double>* itkNotUsed(m2),
384  std::vector<double>* itkNotUsed(m3))
385 {
386  return 0;
387 }
388 
392 template <class TInputImage, class TOutputImage, class TOutputImageDirection, class InterpolatorType>
394 {
395  Superclass::PrintSelf(os, indent);
396  os << indent << "Length: " << m_LengthLine << std::endl;
397  os << indent << "Width: " << m_WidthLine << std::endl;
398 }
400 
401 } // end namespace otb
402 
403 #endif
virtual double ComputeMeasure(std::vector< double > *m1, std::vector< double > *m2, std::vector< double > *m3)
void PrintSelf(std::ostream &os, itk::Indent indent) const override
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
constexpr double CONST_PI
Definition: otbMath.h:49
#define otbMsgDevMacro(x)
Definition: otbMacro.h:116