ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
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"
27
28namespace itk
29{
48template <typename TInputImage, typename TOutputImage, typename TInternalPixel = double>
49class ITK_TEMPLATE_EXPORT GradientVectorFlowImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>
50{
51public:
52 ITK_DISALLOW_COPY_AND_MOVE(GradientVectorFlowImageFilter);
53
56
59
63
65 itkNewMacro(Self);
66
68 itkOverrideGetNameOfClassMacro(GradientVectorFlowImageFilter);
69
71 using InputImageType = TInputImage;
72 using OutputImageType = TOutputImage;
73
74 using IndexType = typename TInputImage::IndexType;
75 using SizeType = typename TInputImage::SizeType;
76 using PixelType = typename TInputImage::PixelType;
77 using OutputImagePointer = typename OutputImageType::Pointer;
78 using RegionType = typename OutputImageType::RegionType;
79
84
86 static constexpr unsigned int ImageDimension = TInputImage::ImageDimension;
87 static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
88
89 using InternalPixelType = TInternalPixel;
94
97
99
100 itkSetObjectMacro(LaplacianFilter, LaplacianFilterType);
101 itkGetModifiableObjectMacro(LaplacianFilter, LaplacianFilterType);
102
103 itkSetMacro(TimeStep, double);
104 itkGetConstMacro(TimeStep, double);
105
106 itkSetMacro(NoiseLevel, double);
107 itkGetConstMacro(NoiseLevel, double);
108
109 itkSetMacro(IterationNum, int);
110 itkGetConstMacro(IterationNum, int);
111
114 itkConceptMacro(OutputHasNumericTraitsCheck,
116
117protected:
119 ~GradientVectorFlowImageFilter() override = default;
120 void
121 PrintSelf(std::ostream & os, Indent indent) const override;
122
123 void
124 GenerateData() override;
125
127 void
129
134 void
136
138 void
140
141private:
142 // parameters;
143 double m_TimeStep{}; // the timestep of each
144 // iteration
145 double m_Steps[Superclass::InputImageDimension]{}; // set to be 1 in all
146 // directions in most cases
147 double m_NoiseLevel{}; // the noise level of the
148 // image
149 int m_IterationNum{}; // the iteration number
150
153
155 InternalImagePointer m_BImage{}; // store the "b" value for every pixel
156
157 typename Superclass::InputImagePointer m_CImage{}; // store the $c_i$ value for
158 // every pixel
159};
160} // end namespace itk
161
162#ifndef ITK_MANUAL_INSTANTIATION
163# include "itkGradientVectorFlowImageFilter.hxx"
164#endif
165
166#endif
InternalImagePointer m_InternalImages[Superclass::InputImageDimension]
void PrintSelf(std::ostream &os, Indent indent) const override
ImageRegionConstIterator< InternalImageType > InternalImageConstIterator
double m_Steps[Superclass::InputImageDimension]
typename LaplacianFilterType::Pointer LaplacianFilterPointer
ImageRegionIterator< OutputImageType > OutputImageIterator
typename OutputImageType::RegionType RegionType
typename OutputImageType::Pointer OutputImagePointer
~GradientVectorFlowImageFilter() override=default
ImageToImageFilter< TInputImage, TOutputImage > Superclass
ImageRegionIterator< InternalImageType > InternalImageIterator
ImageRegionConstIterator< InputImageType > InputImageConstIterator
typename InternalImageType::Pointer InternalImagePointer
ImageRegionIterator< InputImageType > InputImageIterator
LaplacianImageFilter< InternalImageType, InternalImageType > LaplacianFilterType
itk::Image< InternalPixelType, Self::ImageDimension > InternalImageType
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.
typename InputImageType::Pointer InputImagePointer
static constexpr unsigned int InputImageDimension
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.
Implements transparent reference counting.
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....