ITK  4.9.0
Insight Segmentation and Registration Toolkit
itkSatoVesselnessFeatureGenerator.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Insight Segmentation & Registration Toolkit
4  Module: itkSatoVesselnessFeatureGenerator.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 itkSatoVesselnessFeatureGenerator_h
18 #define itkSatoVesselnessFeatureGenerator_h
19 
20 #include "itkFeatureGenerator.h"
21 #include "itkImage.h"
22 #include "itkImageSpatialObject.h"
27 
28 namespace itk
29 {
30 
43 template <unsigned int NDimension>
44 class ITK_EXPORT SatoVesselnessFeatureGenerator : public FeatureGenerator<NDimension>
45 {
46 public:
52 
54  itkNewMacro(Self);
55 
58 
60  itkStaticConstMacro(Dimension, unsigned int, NDimension);
61 
64  typedef signed short InputPixelType;
68  typedef typename Superclass::SpatialObjectType SpatialObjectType;
69 
72  void SetInput( const SpatialObjectType * input );
73  const SpatialObjectType * GetInput() const;
75 
78  const SpatialObjectType * GetFeature() const;
79 
82  itkSetMacro( Sigma, double );
83  itkGetMacro( Sigma, double );
85 
87  itkSetMacro( Alpha1, double );
88  itkGetMacro( Alpha1, double );
90 
92  itkSetMacro( Alpha2, double );
93  itkGetMacro( Alpha2, double );
95 
97  itkSetMacro( UseVesselEnhancingDiffusion, bool );
98  itkGetMacro( UseVesselEnhancingDiffusion, bool );
99  itkBooleanMacro( UseVesselEnhancingDiffusion );
101 
102 protected:
105  void PrintSelf(std::ostream& os, Indent indent) const;
106 
109  void GenerateData ();
110 
111 private:
112  SatoVesselnessFeatureGenerator(const Self&); //purposely not implemented
113  void operator=(const Self&); //purposely not implemented
114 
115  typedef float InternalPixelType;
117 
120 
122 
126 
130 
131  double m_Sigma;
132  double m_Alpha1;
133  double m_Alpha2;
135 };
136 
137 } // end namespace itk
138 
139 #ifndef ITK_MANUAL_INSTANTIATION
140 # include "itkSatoVesselnessFeatureGenerator.hxx"
141 #endif
142 
143 #endif
Light weight base class for most itk classes.
Image< InputPixelType, Dimension > InputImageType
ImageSpatialObject< NDimension, OutputPixelType > OutputImageSpatialObjectType
Class that generates features (typically images) used as input for a segmentation method...
Implementation of the composite pattern.
virtual void SetInput(const DataObjectIdentifierType &key, DataObject *input)
Protected method for setting indexed and named inputs.
Computes the Hessian matrix of an image by convolution with the Second and Cross derivatives of a Gau...
HessianRecursiveGaussianImageFilter< InputImageType > HessianFilterType
Line filter to provide a vesselness measure for tubular objects from the hessian matrix.
ImageSpatialObject< NDimension, InputPixelType > InputImageSpatialObjectType
Generates a feature image by computing the Sato Vesselness measure of the input image.
Implementation of an image as spatial object.
InputImageSpatialObjectType::Pointer InputImageSpatialObjectPointer
Hessian3DToVesselnessMeasureImageFilter< InternalPixelType > VesselnessMeasureFilterType
Control indentation during Print() invocation.
Definition: itkIndent.h:49
VesselnessMeasureFilterType::Pointer m_VesselnessFilter
VesselEnhancingDiffusion3DImageFilter< InputPixelType, Dimension > VesselEnhancingDiffusionFilterType
VesselEnhancingDiffusionFilterType::Pointer m_VesselEnhancingDiffusionFilter
Templated n-dimensional image class.
Definition: itkImage.h:75
Image< InternalPixelType, Dimension > InternalImageType