ITK  4.8.0
Insight Segmentation and Registration Toolkit
itkCannyEdgesDistanceAdvectionFieldFeatureGenerator.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Insight Segmentation & Registration Toolkit
4  Module: itkCannyEdgesDistanceAdvectionFieldFeatureGenerator.h
5  Language: C++
6  Date: $Date$
7  Version: $Revision$
8 
9  Copyright (c) Insight Software Consortium. All rights reserved.
10  See ITKCopyright.txt or http://www.itk.org/HTML/Copyright.htm for details.
11 
12  This software is distributed WITHOUT ANY WARRANTY; without even
13  the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
14  PURPOSE. See the above copyright notices for more information.
15 
16 =========================================================================*/
17 #ifndef __itkCannyEdgesDistanceAdvectionFieldFeatureGenerator_h
18 #define __itkCannyEdgesDistanceAdvectionFieldFeatureGenerator_h
19 
20 #include "itkFeatureGenerator.h"
21 #include "itkImage.h"
22 #include "itkImageSpatialObject.h"
23 #include "itkCastImageFilter.h"
26 #include "itkGradientImageFilter.h"
27 #include "itkMultiplyImageFilter.h"
28 
29 namespace itk
30 {
31 
66 template <unsigned int NDimension>
68 {
69 public:
75 
77  itkNewMacro(Self);
78 
81 
83  itkStaticConstMacro(Dimension, unsigned int, NDimension);
84 
87  typedef signed short InputPixelType;
91  typedef typename Superclass::SpatialObjectType SpatialObjectType;
92 
95  void SetInput( const SpatialObjectType * input );
96  const SpatialObjectType * GetInput() const;
98 
101  const SpatialObjectType * GetFeature() const;
102 
103  itkSetMacro( Sigma, double );
104  itkGetMacro( Sigma, double );
105  itkSetMacro( UpperThreshold, double );
106  itkGetMacro( UpperThreshold, double );
107  itkSetMacro( LowerThreshold, double );
108  itkGetMacro( LowerThreshold, double );
109 
110 protected:
113  void PrintSelf(std::ostream& os, Indent indent) const;
114 
117  void GenerateData ();
118 
119 private:
120  CannyEdgesDistanceAdvectionFieldFeatureGenerator(const Self&); //purposely not implemented
121  void operator=(const Self&); //purposely not implemented
122 
123  typedef float InternalPixelType;
125 
126  typedef CastImageFilter<
132 
136 
141 
142  typedef typename CovariantVectorImageType::PixelType OutputPixelType;
145 
146  typedef MultiplyImageFilter<
150 
156 
159  double m_Sigma;
160 
161 };
162 
163 } // end namespace itk
164 
165 #ifndef ITK_MANUAL_INSTANTIATION
166 # include "itkCannyEdgesDistanceAdvectionFieldFeatureGenerator.hxx"
167 #endif
168 
169 #endif
Light weight base class for most itk classes.
MultiplyImageFilter< CovariantVectorImageType, InternalImageType, CovariantVectorImageType > MultiplyFilterType
SignedMaurerDistanceMapImageFilter< InternalImageType, InternalImageType > DistanceMapFilterType
This filter calculates the Euclidean distance transform of a binary image in linear time for arbitrar...
Class that generates features (typically images) used as input for a segmentation method...
virtual void SetInput(const DataObjectIdentifierType &key, DataObject *input)
Protected method for setting indexed and named inputs.
Implementation of an image as spatial object.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Pixel-wise multiplication of two images.
GradientImageFilter< InternalImageType, InternalPixelType, InternalPixelType > GradientFilterType
Generates an advection feature field by computing the distance map to the canny edges in the image an...
CannyEdgeDetectionRecursiveGaussianImageFilter< InternalImageType, InternalImageType > CannyEdgeFilterType
Templated n-dimensional image class.
Definition: itkImage.h:75
Casts input pixels to output pixel type.
Computes the gradient of an image using directional derivatives.