OTB  10.0.0
Orfeo Toolbox
otbPixelSuppressionByDirectionImageFilter.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 otbPixelSuppressionByDirectionImageFilter_hxx
22 #define otbPixelSuppressionByDirectionImageFilter_hxx
23 
24 #include "otbMacro.h" //for
26 
27 #include "itkDataObject.h"
28 #include "itkConstNeighborhoodIterator.h"
29 #include "itkNeighborhoodInnerProduct.h"
30 #include "itkImageRegionIterator.h"
31 #include "itkNeighborhoodAlgorithm.h"
32 #include "itkConstantBoundaryCondition.h"
33 #include "itkOffset.h"
34 #include "itkProgressReporter.h"
35 #include "otbMath.h"
36 
37 namespace otb
38 {
39 
43 template <class TInputImage, class TOutputImage>
45 {
46  m_Radius.Fill(1);
47  m_AngularBeam = static_cast<double>(0.);
48  this->DynamicMultiThreadingOn();
49 }
51 
52 template <class TInputImage, class TOutputImage>
54 {
55  this->SetInput(0, input);
56 }
57 
58 template <class TInputImage, class TOutputImage>
60 {
61  this->SetInput(1, input);
62 }
63 
64 template <class TInputImage, class TOutputImage>
67 {
68  if (this->GetNumberOfInputs() < 1)
69  {
70  return nullptr;
71  }
72 
73  return static_cast<const TInputImage*>(this->GetInput(0));
74 }
75 
76 template <class TInputImage, class TOutputImage>
79 {
80  if (this->GetNumberOfInputs() < 1)
81  {
82  return nullptr;
83  }
84 
85  return static_cast<const TInputImage*>(this->GetInput(1));
86 }
87 
88 template <class TInputImage, class TOutputImage>
90 {
91  // call the superclass' implementation of this method
92  Superclass::GenerateInputRequestedRegion();
93 
94  // get pointers to the input and output
95  typename Superclass::InputImagePointer inputPtr = const_cast<TInputImage*>(this->GetInputImageDirection());
96  typename Superclass::OutputImagePointer outputPtr = this->GetOutput();
97 
98  if (!inputPtr || !outputPtr)
99  {
100  return;
101  }
102 
103  // get a copy of the input requested region (should equal the output
104  // requested region)
105  typename TInputImage::RegionType inputRequestedRegion;
106  inputRequestedRegion = inputPtr->GetRequestedRegion();
107 
108  // pad the input requested region by the operator radius
109  inputRequestedRegion.PadByRadius(m_Radius);
110 
111  // crop the input requested region at the input's largest possible region
112  if (inputRequestedRegion.Crop(inputPtr->GetLargestPossibleRegion()))
113  {
114  inputPtr->SetRequestedRegion(inputRequestedRegion);
115  return;
116  }
117  else
118  {
119  // Couldn't crop the region (requested region is outside the largest
120  // possible region). Throw an exception.
121 
122  // store what we tried to request (prior to trying to crop)
123  inputPtr->SetRequestedRegion(inputRequestedRegion);
124 
125  // build an exception
126  itk::InvalidRequestedRegionError e(__FILE__, __LINE__);
127  std::ostringstream msg;
128  msg << static_cast<const char*>(this->GetNameOfClass()) << "::GenerateInputRequestedRegion()";
129  e.SetLocation(msg.str());
130  e.SetDescription("Requested region is (at least partially) outside the largest possible region.");
131  e.SetDataObject(inputPtr);
132  throw e;
133  }
134 }
135 
136 template <class TInputImage, class TOutputImage>
138 {
139 
140  itk::ConstantBoundaryCondition<InputImageType> cbc;
141  const InputPixelType cvalue = 255;
142  cbc.SetConstant(cvalue);
143 
144  itk::ConstNeighborhoodIterator<InputImageType> bit;
145  itk::ImageRegionConstIterator<InputImageType> itin;
146  itk::ImageRegionIterator<OutputImageType> itout;
147 
148  // Allocate output
149  typename OutputImageType::Pointer output = this->GetOutput();
150  typename InputImageType::ConstPointer input = this->GetInputImage();
151  typename InputImageType::ConstPointer inputDirection = this->GetInputImageDirection();
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(inputDirection, outputRegionForThread, m_Radius);
159 
160  // typename TInputImage::IndexType bitIndex;
161 
162  //---------------------------------------------------------------------------
163 
164  InputPixelType PixelValue;
165 
166  // Location of the central pixel in the input image
167  // int Xc, Yc;
168 
169  // Pixel location in the system axis of the region
170  int x, y;
171 
172  // Pixel location in the system axis of the region after rotation of theta
173  // where theta is the direction of the cantral pixel
174 
175  // Pixel Direction
176  double ThetaXcYc, Thetaxtyt;
177 
178  //---------------------------------------------------------------------------
179 
180  // Process each of the boundary faces. These are N-d regions which border
181  // the edge of the buffer.
182  for (fit = faceList.begin(); fit != faceList.end(); ++fit)
183  {
184  bit = itk::ConstNeighborhoodIterator<InputImageType>(m_Radius, inputDirection, *fit);
185 
186  itin = itk::ImageRegionConstIterator<InputImageType>(input, *fit);
187  itout = itk::ImageRegionIterator<OutputImageType>(output, *fit);
188 
189  bit.OverrideBoundaryCondition(&cbc);
190  bit.GoToBegin();
191 
192  while (!bit.IsAtEnd())
193  {
194 
195  /*// Location of the central pixel of the region in the input image
196  bitIndex = bit.GetIndex();
197 
198  Xc = bitIndex[0];
199  Yc = bitIndex[1];
200  */
201  // Get Pixel Direction from the image of directions
202  ThetaXcYc = static_cast<double>(bit.GetCenterPixel());
203 
204  // Pixel intensity in the input image
205  PixelValue = itin.Get();
206 
207  bool IsLine = false;
208 
209  typename itk::ConstNeighborhoodIterator<InputImageType>::OffsetType off;
210  // Loop on the region
211  for (unsigned int i = 0; i < 2 * m_Radius[0] + 1; ++i)
212  for (unsigned int j = 0; j < 2 * m_Radius[1] + 1; ++j)
213  {
214 
215  off[0] = i - m_Radius[0];
216  off[1] = j - m_Radius[1];
217 
218  x = off[0];
219  y = off[1];
220 
221  // No calculation on the central pixel
222  if ((x == 0) && (y == 0))
223  continue;
224 
225  Thetaxtyt = std::atan2(static_cast<double>(y), static_cast<double>(x)); // result is [-PI, PI]
226  while (Thetaxtyt < 0)
227  Thetaxtyt = CONST_PI + Thetaxtyt; // Theta is now [0, PI] as is
228  // the result of detectors
229  while (Thetaxtyt > CONST_PI_2)
230  Thetaxtyt = Thetaxtyt - CONST_PI; // Theta is now [-PI/2, PI/2]
231 
232  if ((std::abs(std::cos(Thetaxtyt - ThetaXcYc)) >= std::cos(m_AngularBeam)) // this
233  // pixel
234  // is
235  // in
236  // the
237  // angular beam
238  && (std::abs(std::cos(bit.GetPixel(off) - ThetaXcYc)) >= std::cos(m_AngularBeam))) // and
239  // its
240  // direction
241  // is
242  // also
243  // in
244  // the beam
245  {
246  IsLine = true;
247  continue;
248  }
249  }
250 
251  // end of the loop on the pixels of the region
252 
253  // Assignment of this value to the output pixel
254  if (IsLine == true)
255  {
256  itout.Set(static_cast<OutputPixelType>(PixelValue));
257  }
258  else
259  {
260  itout.Set(static_cast<OutputPixelType>(0.));
261  }
262 
263  ++bit;
264  ++itin;
265  ++itout;
266  }
267  }
268 }
269 
273 template <class TInputImage, class TOutput>
274 void PixelSuppressionByDirectionImageFilter<TInputImage, TOutput>::PrintSelf(std::ostream& os, itk::Indent indent) const
275 {
276  Superclass::PrintSelf(os, indent);
277  os << indent << "Radius: " << m_Radius << std::endl;
278 }
280 
281 } // end namespace otb
282 
283 #endif
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_2
Definition: otbMath.h:50
constexpr double CONST_PI
Definition: otbMath.h:49