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