ITK  5.2.0
Insight Toolkit
itkGradientVectorFlowImageFilter.h
Go to the documentation of this file.
1 /*=========================================================================
2  *
3  * Copyright NumFOCUS
4  *
5  * Licensed under the Apache License, Version 2.0 (the "License");
6  * you may not use this file except in compliance with the License.
7  * You may obtain a copy of the License at
8  *
9  * http://www.apache.org/licenses/LICENSE-2.0.txt
10  *
11  * Unless required by applicable law or agreed to in writing, software
12  * distributed under the License is distributed on an "AS IS" BASIS,
13  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
14  * See the License for the specific language governing permissions and
15  * limitations under the License.
16  *
17  *=========================================================================*/
18 #ifndef itkGradientVectorFlowImageFilter_h
19 #define itkGradientVectorFlowImageFilter_h
20 
21 #include "vnl/vnl_matrix_fixed.h"
22 #include "itkMath.h"
23 #include "itkImage.h"
24 #include "itkVector.h"
26 #include "itkImageRegionIterator.h"
28 
29 namespace itk
30 {
49 template <typename TInputImage, typename TOutputImage, typename TInternalPixel = double>
50 class ITK_TEMPLATE_EXPORT GradientVectorFlowImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
51 {
52 public:
53  ITK_DISALLOW_COPY_AND_MOVE(GradientVectorFlowImageFilter);
54 
57 
60 
64 
66  itkNewMacro(Self);
67 
70 
72  using InputImageType = TInputImage;
73  using OutputImageType = TOutputImage;
74 
76  using SizeType = typename TInputImage::SizeType;
77  using PixelType = typename TInputImage::PixelType;
78  using OutputImagePointer = typename OutputImageType::Pointer;
80 
85 
87  static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
88  static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
89 
90  using InternalPixelType = TInternalPixel;
95 
98 
101  itkSetObjectMacro(LaplacianFilter, LaplacianFilterType);
102  itkGetModifiableObjectMacro(LaplacianFilter, LaplacianFilterType);
103 
104  itkSetMacro(TimeStep, double);
105  itkGetConstMacro(TimeStep, double);
106 
107  itkSetMacro(NoiseLevel, double);
108  itkGetConstMacro(NoiseLevel, double);
109 
110  itkSetMacro(IterationNum, int);
111  itkGetConstMacro(IterationNum, int);
112 
113 #ifdef ITK_USE_CONCEPT_CHECKING
114  // Begin concept checking
117  itkConceptMacro(OutputHasNumericTraitsCheck,
119  // End concept checking
120 #endif
121 
122 protected:
124  ~GradientVectorFlowImageFilter() override = default;
125  void
126  PrintSelf(std::ostream & os, Indent indent) const override;
127 
128  void
129  GenerateData() override;
130 
132  void
133  InitInterImage();
134 
139  void
140  UpdateInterImage();
141 
143  void
144  UpdatePixels();
145 
146 private:
147  // parameters;
148  double m_TimeStep; // the timestep of each
149  // iteration
150  double m_Steps[Superclass::InputImageDimension]; // set to be 1 in all
151  // directions in most cases
152  double m_NoiseLevel; // the noise level of the
153  // image
154  int m_IterationNum; // the iteration number
155 
158 
159  InternalImagePointer m_InternalImages[Superclass::InputImageDimension];
160  InternalImagePointer m_BImage; // store the "b" value for every pixel
161 
162  typename Superclass::InputImagePointer m_CImage; // store the $c_i$ value for
163  // every pixel
164 };
165 } // end namespace itk
166 
167 #ifndef ITK_MANUAL_INSTANTIATION
168 # include "itkGradientVectorFlowImageFilter.hxx"
169 #endif
170 
171 #endif
itk::GradientVectorFlowImageFilter::m_LaplacianFilter
LaplacianFilterPointer m_LaplacianFilter
Definition: itkGradientVectorFlowImageFilter.h:156
itk::Concept::HasNumericTraits
Definition: itkConceptChecking.h:714
itk::ImageSource::OutputImagePointer
typename OutputImageType::Pointer OutputImagePointer
Definition: itkImageSource.h:91
itkImageRegionConstIteratorWithIndex.h
itk::GradientVectorFlowImageFilter::SizeType
typename TInputImage::SizeType SizeType
Definition: itkGradientVectorFlowImageFilter.h:76
itk::GradientVectorFlowImageFilter::IndexType
typename TInputImage::IndexType IndexType
Definition: itkGradientVectorFlowImageFilter.h:75
itk::LaplacianImageFilter
This filter computes the Laplacian of a scalar-valued image.
Definition: itkLaplacianImageFilter.h:63
itk::GradientVectorFlowImageFilter::PixelType
typename TInputImage::PixelType PixelType
Definition: itkGradientVectorFlowImageFilter.h:77
itk::GradientVectorFlowImageFilter::m_BImage
InternalImagePointer m_BImage
Definition: itkGradientVectorFlowImageFilter.h:160
itk::GTest::TypedefsAndConstructors::Dimension2::SizeType
ImageBaseType::SizeType SizeType
Definition: itkGTestTypedefsAndConstructors.h:49
itk::GradientVectorFlowImageFilter
This class computes a diffusion of the gradient vectors for graylevel or binary edge map derive from ...
Definition: itkGradientVectorFlowImageFilter.h:50
itkImage.h
itk::SmartPointer< Self >
itkImageRegionIterator.h
itk::Indent
Control indentation during Print() invocation.
Definition: itkIndent.h:49
itk::GradientVectorFlowImageFilter::m_IntermediateImage
Superclass::InputImagePointer m_IntermediateImage
Definition: itkGradientVectorFlowImageFilter.h:157
itkLaplacianImageFilter.h
itk::Concept::SameDimension
Definition: itkConceptChecking.h:694
itk::ImageRegionIterator
A multi-dimensional iterator templated over image type that walks a region of pixels.
Definition: itkImageRegionIterator.h:80
itk::GTest::TypedefsAndConstructors::Dimension2::IndexType
ImageBaseType::IndexType IndexType
Definition: itkGTestTypedefsAndConstructors.h:50
itk::ImageToImageFilter
Base class for filters that take an image as input and produce an image as output.
Definition: itkImageToImageFilter.h:108
itk::ImageSource
Base class for all process objects that output image data.
Definition: itkImageSource.h:67
itk::GradientVectorFlowImageFilter::m_IterationNum
int m_IterationNum
Definition: itkGradientVectorFlowImageFilter.h:154
itk::GradientVectorFlowImageFilter::InternalPixelType
TInternalPixel InternalPixelType
Definition: itkGradientVectorFlowImageFilter.h:90
itk::GradientVectorFlowImageFilter::RegionType
typename OutputImageType::RegionType RegionType
Definition: itkGradientVectorFlowImageFilter.h:79
itk::ImageToImageFilter::InputImagePointer
typename InputImageType::Pointer InputImagePointer
Definition: itkImageToImageFilter.h:130
itk::GTest::TypedefsAndConstructors::Dimension2::RegionType
ImageBaseType::RegionType RegionType
Definition: itkGTestTypedefsAndConstructors.h:54
itk::GradientVectorFlowImageFilter::m_CImage
Superclass::InputImagePointer m_CImage
Definition: itkGradientVectorFlowImageFilter.h:162
itk::ImageToImageFilter::InputImageType
TInputImage InputImageType
Definition: itkImageToImageFilter.h:129
itk::GradientVectorFlowImageFilter::LaplacianFilterPointer
typename LaplacianFilterType::Pointer LaplacianFilterPointer
Definition: itkGradientVectorFlowImageFilter.h:97
itk::GradientVectorFlowImageFilter::m_TimeStep
double m_TimeStep
Definition: itkGradientVectorFlowImageFilter.h:148
itk::GradientVectorFlowImageFilter::m_NoiseLevel
double m_NoiseLevel
Definition: itkGradientVectorFlowImageFilter.h:152
itkConceptMacro
#define itkConceptMacro(name, concept)
Definition: itkConceptChecking.h:65
itk
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
Definition: itkAnnulusOperator.h:24
itk::ProcessObject
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
Definition: itkProcessObject.h:138
itkVector.h
itk::ImageRegionConstIterator
A multi-dimensional iterator templated over image type that walks a region of pixels.
Definition: itkImageRegionConstIterator.h:109
itk::Image
Templated n-dimensional image class.
Definition: itkImage.h:86
itk::GradientVectorFlowImageFilter::InternalImagePointer
typename InternalImageType::Pointer InternalImagePointer
Definition: itkGradientVectorFlowImageFilter.h:92
itkMath.h
itk::ImageSource::OutputImageType
TOutputImage OutputImageType
Definition: itkImageSource.h:90