ITK  5.0.0
Insight Segmentation and Registration Toolkit
itkFrangiTubularnessFeatureGenerator.h
Go to the documentation of this file.
1 /*=========================================================================
2 
3  Program: Insight Segmentation & Registration Toolkit
4  Module: itkFrangiTubularnessFeatureGenerator.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 itkFrangiTubularnessFeatureGenerator_h
18 #define itkFrangiTubularnessFeatureGenerator_h
19 
20 #include "itkFeatureGenerator.h"
21 #include "itkImage.h"
22 #include "itkImageSpatialObject.h"
27 
28 namespace itk
29 {
30 
41 template <unsigned int NDimension>
42 class ITK_TEMPLATE_EXPORT FrangiTubularnessFeatureGenerator : public FeatureGenerator<NDimension>
43 {
44 public:
45  ITK_DISALLOW_COPY_AND_ASSIGN(FrangiTubularnessFeatureGenerator);
46 
52 
54  itkNewMacro(Self);
55 
58 
60  static constexpr unsigned int Dimension = NDimension;
61 
64  using InputPixelType = signed short;
68  using SpatialObjectType = typename Superclass::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( SheetnessNormalization, double );
88  itkGetMacro( SheetnessNormalization, double );
90 
92  itkSetMacro( BloobinessNormalization, double );
93  itkGetMacro( BloobinessNormalization, double );
95 
97  itkSetMacro( NoiseNormalization, double );
98  itkGetMacro( NoiseNormalization, double );
100 
101 protected:
104  void PrintSelf(std::ostream& os, Indent indent) const override;
105 
108  void GenerateData () override;
109 
110 private:
111  using InternalPixelType = float;
113 
116 
118 
121  using HessianPixelType = typename HessianImageType::PixelType;
122 
125 
127 
129 
133 
134  double m_Sigma;
138 };
139 
140 } // end namespace itk
141 
142 #ifndef ITK_MANUAL_INSTANTIATION
143 # include "itkFrangiTubularnessFeatureGenerator.hxx"
144 #endif
145 
146 #endif
Light weight base class for most itk classes.
Computes the eigen-values of every input symmetric matrix pixel.
typename HessianFilterType::OutputImageType HessianImageType
Simulate a standard C array with copy semnatics.
Definition: itkFixedArray.h:51
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.
Generates a feature image by computing measures based on the Hessian Eigenvalues. ...
Computes the Hessian matrix of an image by convolution with the Second and Cross derivatives of a Gau...
Computes a measure of CrestLines from the Hessian Eigenvalues.
typename InputImageSpatialObjectType::Pointer InputImageSpatialObjectPointer
Implementation of an image as spatial object.
Control indentation during Print() invocation.
Definition: itkIndent.h:49
Templated n-dimensional image class.
Definition: itkImage.h:75