21 #ifndef otbTouziEdgeDetectorImageFilter_hxx
22 #define otbTouziEdgeDetectorImageFilter_hxx
26 #include "itkDataObject.h"
27 #include "itkConstNeighborhoodIterator.h"
28 #include "itkNeighborhoodInnerProduct.h"
29 #include "itkImageRegionIterator.h"
30 #include "itkNeighborhoodAlgorithm.h"
31 #include "itkProgressReporter.h"
41 template <
class TInputImage,
class TOutputImage,
class TOutputImageDirection>
45 this->DynamicMultiThreadingOn();
49 template <
class TInputImage,
class TOutputImage,
class TOutputImageDirection>
53 Superclass::GenerateInputRequestedRegion();
59 if (!inputPtr || !outputPtr)
66 typename TInputImage::RegionType inputRequestedRegion;
67 inputRequestedRegion = inputPtr->GetRequestedRegion();
70 inputRequestedRegion.PadByRadius(m_Radius);
73 if (inputRequestedRegion.Crop(inputPtr->GetLargestPossibleRegion()))
75 inputPtr->SetRequestedRegion(inputRequestedRegion);
84 inputPtr->SetRequestedRegion(inputRequestedRegion);
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);
102 template <
class TInputImage,
class TOutputImage,
class TOutputImageDirection>
106 typename OutputImageDirectionType::RegionType region;
107 typename OutputImageType::Pointer output = this->GetOutput();
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();
119 template <
class TInputImage,
class TOutputImage,
class TOutputImageDirection>
123 itk::ZeroFluxNeumannBoundaryCondition<InputImageType> nbc;
124 itk::ConstNeighborhoodIterator<InputImageType> bit;
125 itk::ImageRegionIterator<OutputImageType> it;
126 itk::ImageRegionIterator<OutputImageType> it_dir;
129 typename OutputImageType::Pointer output = this->GetOutput();
130 typename InputImageType::ConstPointer input = this->GetInput();
131 typename OutputImageDirectionType::Pointer outputDir = this->GetOutputDirection();
134 typename itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<InputImageType>::FaceListType faceList;
135 typename itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<InputImageType>::FaceListType::iterator fit;
137 itk::NeighborhoodAlgorithm::ImageBoundaryFacesCalculator<InputImageType> bC;
138 faceList = bC(input, outputRegionForThread, m_Radius);
141 typename TInputImage::IndexType bitIndex;
146 const int NB_DIR = 4;
148 const int NB_REGION = 2;
150 double Theta[NB_DIR];
158 double Sum[NB_DIR][NB_REGION];
164 double R_theta[NB_DIR];
165 double Sum_R_theta = 0.;
169 double Dir_contour = 0.;
184 for (fit = faceList.begin(); fit != faceList.end(); ++fit)
189 bit = itk::ConstNeighborhoodIterator<InputImageType>(m_Radius, input, *fit);
190 unsigned int neighborhoodSize = bit.Size();
192 it = itk::ImageRegionIterator<OutputImageType>(output, *fit);
193 it_dir = itk::ImageRegionIterator<OutputImageDirectionType>(outputDir, *fit);
195 bit.OverrideBoundaryCondition(&nbc);
198 while (!bit.IsAtEnd())
202 bitIndex = bit.GetIndex();
208 for (
int dir = 0; dir < NB_DIR; ++dir)
210 for (
int m = 0; m < NB_REGION; m++)
219 for (i = 0; i < neighborhoodSize; ++i)
222 bitIndex = bit.GetIndex(i);
230 Sum[0][0] +=
static_cast<double>(bit.GetPixel(i));
232 Sum[0][1] +=
static_cast<double>(bit.GetPixel(i));
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));
242 Sum[2][0] +=
static_cast<double>(bit.GetPixel(i));
244 Sum[2][1] +=
static_cast<double>(bit.GetPixel(i));
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));
255 for (
int dir = 0; dir < NB_DIR; ++dir)
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));
262 if ((M1 != 0) && (M2 != 0))
263 R_theta[dir] =
static_cast<double>(1 - std::min((M1 / M2), (M2 / M1)));
268 R_contour =
static_cast<double>(std::max(R_contour, R_theta[dir]));
276 Dir_contour += sign * Theta[dir] * R_theta[dir];
277 Sum_R_theta += R_theta[dir];
285 if (Sum_R_theta != 0.)
286 Dir_contour = Dir_contour / Sum_R_theta;
301 template <
class TInputImage,
class TOutputImage,
class TOutputImageDirection>
304 Superclass::PrintSelf(os, indent);
305 os << indent <<
"Radius: " << m_Radius << std::endl;
OutputImageType::Pointer OutputImagePointer
InputImageType::Pointer InputImagePointer
TouziEdgeDetectorImageFilter()
void GenerateInputRequestedRegion() override
OutputImageType::RegionType OutputImageRegionType
void BeforeThreadedGenerateData() override
void PrintSelf(std::ostream &os, itk::Indent indent) const override
void DynamicThreadedGenerateData(const OutputImageRegionType &outputRegionForThread) override
OutputImageDirectionType::PixelType OutputPixelDirectionType
OutputImageType::PixelType OutputPixelType
Superclass::OutputImageDirectionType OutputImageDirectionType
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
constexpr double CONST_PI_4
constexpr double CONST_PI_2
constexpr double CONST_PI