ITK  6.0.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 * https://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"
28
29namespace itk
30{
49template <typename TInputImage, typename TOutputImage, typename TInternalPixel = double>
50class ITK_TEMPLATE_EXPORT GradientVectorFlowImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
51{
52public:
53 ITK_DISALLOW_COPY_AND_MOVE(GradientVectorFlowImageFilter);
54
57
60
64
66 itkNewMacro(Self);
67
69 itkOverrideGetNameOfClassMacro(GradientVectorFlowImageFilter);
70
72 using InputImageType = TInputImage;
73 using OutputImageType = TOutputImage;
74
77 using PixelType = typename TInputImage::PixelType;
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
122protected:
124 ~GradientVectorFlowImageFilter() override = default;
125 void
126 PrintSelf(std::ostream & os, Indent indent) const override;
127
128 void
129 GenerateData() override;
130
132 void
134
139 void
141
143 void
145
146private:
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
156 LaplacianFilterPointer m_LaplacianFilter{};
157 typename Superclass::InputImagePointer m_IntermediateImage{};
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
This class computes a diffusion of the gradient vectors for graylevel or binary edge map derive from ...
void PrintSelf(std::ostream &os, Indent indent) const override
typename LaplacianFilterType::Pointer LaplacianFilterPointer
typename OutputImageType::RegionType RegionType
~GradientVectorFlowImageFilter() override=default
typename InternalImageType::Pointer InternalImagePointer
A multi-dimensional iterator templated over image type that walks a region of pixels.
A multi-dimensional iterator templated over image type that walks a region of pixels.
Base class for all process objects that output image data.
TOutputImage OutputImageType
typename OutputImageType::Pointer OutputImagePointer
Base class for filters that take an image as input and produce an image as output.
typename InputImageType::Pointer InputImagePointer
Templated n-dimensional image class.
Definition: itkImage.h:89
Control indentation during Print() invocation.
Definition: itkIndent.h:50
This filter computes the Laplacian of a scalar-valued image.
The base class for all process objects (source, filters, mappers) in the Insight data processing pipe...
SmartPointer< Self > Pointer
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....