OTB  10.0.0
Orfeo Toolbox
otbTouziEdgeDetectorImageFilter.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 otbTouziEdgeDetectorImageFilter_hxx
22 #define otbTouziEdgeDetectorImageFilter_hxx
23 
25 
26 #include "itkDataObject.h"
27 #include "itkConstNeighborhoodIterator.h"
28 #include "itkNeighborhoodInnerProduct.h"
29 #include "itkImageRegionIterator.h"
30 #include "itkNeighborhoodAlgorithm.h"
31 #include "itkProgressReporter.h"
32 #include "otbMacro.h"
33 #include "otbMath.h"
34 
35 namespace otb
36 {
37 
41 template <class TInputImage, class TOutputImage, class TOutputImageDirection>
43 {
44  m_Radius.Fill(1);
45  this->DynamicMultiThreadingOn();
46 }
48 
49 template <class TInputImage, class TOutputImage, class TOutputImageDirection>
51 {
52  // call the superclass' implementation of this method
53  Superclass::GenerateInputRequestedRegion();
54 
55  // get pointers to the input and output
56  typename Superclass::InputImagePointer inputPtr = const_cast<TInputImage*>(this->GetInput());
57  typename Superclass::OutputImagePointer outputPtr = this->GetOutput();
58 
59  if (!inputPtr || !outputPtr)
60  {
61  return;
62  }
63 
64  // get a copy of the input requested region (should equal the output
65  // requested region)
66  typename TInputImage::RegionType inputRequestedRegion;
67  inputRequestedRegion = inputPtr->GetRequestedRegion();
68 
69  // pad the input requested region by the operator radius
70  inputRequestedRegion.PadByRadius(m_Radius);
71 
72  // crop the input requested region at the input's largest possible region
73  if (inputRequestedRegion.Crop(inputPtr->GetLargestPossibleRegion()))
74  {
75  inputPtr->SetRequestedRegion(inputRequestedRegion);
76  return;
77  }
78  else
79  {
80  // Couldn't crop the region (requested region is outside the largest
81  // possible region). Throw an exception.
82 
83  // store what we tried to request (prior to trying to crop)
84  inputPtr->SetRequestedRegion(inputRequestedRegion);
85 
86  // build an exception
87  itk::InvalidRequestedRegionError e(__FILE__, __LINE__);
88  std::ostringstream msg;
89  msg << static_cast<const char*>(this->GetNameOfClass()) << "::GenerateInputRequestedRegion()";
90  e.SetLocation(msg.str());
91  e.SetDescription("Requested region is (at least partially) outside the largest possible region.");
92  e.SetDataObject(inputPtr);
93  throw e;
94  }
95 }
96 
102 template <class TInputImage, class TOutputImage, class TOutputImageDirection>
104 {
105 
106  typename OutputImageDirectionType::RegionType region;
107  typename OutputImageType::Pointer output = this->GetOutput();
108 
109  OutputImageDirectionType* direction = this->GetOutputDirection();
110 
111  region.SetSize(output->GetRequestedRegion().GetSize());
112  region.SetIndex(output->GetRequestedRegion().GetIndex());
113  direction->SetRegions(region);
114  direction->SetOrigin(output->GetOrigin());
115  direction->SetSignedSpacing(output->GetSignedSpacing());
116  direction->Allocate();
117 }
118 
119 template <class TInputImage, class TOutputImage, class TOutputImageDirection>
121 {
122  unsigned int i;
123  itk::ZeroFluxNeumannBoundaryCondition<InputImageType> nbc;
124  itk::ConstNeighborhoodIterator<InputImageType> bit;
125  itk::ImageRegionIterator<OutputImageType> it;
126  itk::ImageRegionIterator<OutputImageType> it_dir;
127 
128  // Allocate output
129  typename OutputImageType::Pointer output = this->GetOutput();
130  typename InputImageType::ConstPointer input = this->GetInput();
131  typename OutputImageDirectionType::Pointer outputDir = this->GetOutputDirection();
132 
133  // Find the data-set boundary "faces"
134  typename itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<InputImageType>::FaceListType faceList;
135  typename itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<InputImageType>::FaceListType::iterator fit;
136 
137  itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<InputImageType> bC;
138  faceList = bC(input, outputRegionForThread, m_Radius);
139 
140 
141  typename TInputImage::IndexType bitIndex;
142 
143  // Initializations
144  // ---------------
145  // Number of direction
146  const int NB_DIR = 4;
147  // Number of region of the filter
148  const int NB_REGION = 2;
149  // Definition of the 4 directions
150  double Theta[NB_DIR];
151 
152  Theta[0] = 0.;
153  Theta[1] = CONST_PI_4;
154  Theta[2] = CONST_PI_2;
155  Theta[3] = 3 * CONST_PI / 4.;
156 
157  // contains for the 4 directions the sum of the pixels belonging to each region
158  double Sum[NB_DIR][NB_REGION];
159  // Mean of region 1
160  double M1;
161  // Mean of region 2
162  double M2;
163  // Result of the filter for each direction
164  double R_theta[NB_DIR];
165  double Sum_R_theta = 0.;
166  // Intensity of the contour
167  double R_contour;
168  // Direction of the contour
169  double Dir_contour = 0.;
170  // sign of the contour
171  int sign;
172 
173  // Pixel location in the input image
174  int x;
175  int y;
176  // Location of the central pixel in the input image
177  int xc;
178  int yc;
179 
180  int cpt = 0;
181 
182  // Process each of the boundary faces. These are N-d regions which border
183  // the edge of the buffer.
184  for (fit = faceList.begin(); fit != faceList.end(); ++fit)
185  {
186 
187  cpt += 1;
188 
189  bit = itk::ConstNeighborhoodIterator<InputImageType>(m_Radius, input, *fit);
190  unsigned int neighborhoodSize = bit.Size();
191 
192  it = itk::ImageRegionIterator<OutputImageType>(output, *fit);
193  it_dir = itk::ImageRegionIterator<OutputImageDirectionType>(outputDir, *fit);
194 
195  bit.OverrideBoundaryCondition(&nbc);
196  bit.GoToBegin();
197 
198  while (!bit.IsAtEnd())
199  {
200 
201  // Location of the pixel central
202  bitIndex = bit.GetIndex();
203 
204  xc = bitIndex[0];
205  yc = bitIndex[1];
206 
207  // Initializations
208  for (int dir = 0; dir < NB_DIR; ++dir)
209  {
210  for (int m = 0; m < NB_REGION; m++)
211  Sum[dir][m] = 0.;
212  }
213 
214  R_contour = -1;
215  Dir_contour = 0.;
216  Sum_R_theta = 0.;
217 
218  // Loop on pixels of the filter
219  for (i = 0; i < neighborhoodSize; ++i)
220  {
221 
222  bitIndex = bit.GetIndex(i);
223  x = bitIndex[0];
224  y = bitIndex[1];
225 
226  // We determine for each direction with which region the pixel belongs.
227 
228  // Horizontal direction
229  if (y < yc)
230  Sum[0][0] += static_cast<double>(bit.GetPixel(i));
231  else if (y > yc)
232  Sum[0][1] += static_cast<double>(bit.GetPixel(i));
233 
234  // Diagonal direction 1
235  if ((y - yc) < (x - xc))
236  Sum[1][0] += static_cast<double>(bit.GetPixel(i));
237  else if ((y - yc) > (x - xc))
238  Sum[1][1] += static_cast<double>(bit.GetPixel(i));
239 
240  // Vertical direction
241  if (x > xc)
242  Sum[2][0] += static_cast<double>(bit.GetPixel(i));
243  else if (x < xc)
244  Sum[2][1] += static_cast<double>(bit.GetPixel(i));
245 
246  // Diagonal direction 2
247  if ((y - yc) > -(x - xc))
248  Sum[3][0] += static_cast<double>(bit.GetPixel(i));
249  else if ((y - yc) < -(x - xc))
250  Sum[3][1] += static_cast<double>(bit.GetPixel(i));
251 
252  } // end of the loop on pixels of the filter
253 
254  // Loop on the 4 directions
255  for (int dir = 0; dir < NB_DIR; ++dir)
256  {
257  // Calculation of the mean of the 2 regions
258  M1 = Sum[dir][0] / static_cast<double>(m_Radius[0] * (2 * m_Radius[0] + 1));
259  M2 = Sum[dir][1] / static_cast<double>(m_Radius[0] * (2 * m_Radius[0] + 1));
260 
261  // Calculation of the intensity of the contour
262  if ((M1 != 0) && (M2 != 0))
263  R_theta[dir] = static_cast<double>(1 - std::min((M1 / M2), (M2 / M1)));
264  else
265  R_theta[dir] = 0.;
266 
267  // Determination of the maximum intensity of the contour
268  R_contour = static_cast<double>(std::max(R_contour, R_theta[dir]));
269 
270  // Determination of the sign of contour
271  if (M2 > M1)
272  sign = +1;
273  else
274  sign = -1;
275 
276  Dir_contour += sign * Theta[dir] * R_theta[dir];
277  Sum_R_theta += R_theta[dir];
278 
279  } // end of the loop on the directions
280 
281  // Assignment of this value to the output pixel
282  it.Set(static_cast<OutputPixelType>(R_contour));
283 
284  // Determination of the direction of the contour
285  if (Sum_R_theta != 0.)
286  Dir_contour = Dir_contour / Sum_R_theta;
287 
288  // Assignment of this value to the "outputdir" pixel
289  it_dir.Set(static_cast<OutputPixelDirectionType>(Dir_contour));
290 
291  ++bit;
292  ++it;
293  ++it_dir;
294  }
295  }
296 }
297 
301 template <class TInputImage, class TOutputImage, class TOutputImageDirection>
303 {
304  Superclass::PrintSelf(os, indent);
305  os << indent << "Radius: " << m_Radius << std::endl;
306 }
308 
309 } // end namespace otb
310 
311 #endif
void PrintSelf(std::ostream &os, itk::Indent indent) const override
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
OutputImageDirectionType::PixelType OutputPixelDirectionType
Superclass::OutputImageDirectionType OutputImageDirectionType
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
constexpr double CONST_PI_4
Definition: otbMath.h:51
constexpr double CONST_PI_2
Definition: otbMath.h:50
constexpr double CONST_PI
Definition: otbMath.h:49