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"
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
75 using IndexType = typename TInputImage::IndexType;
76 using SizeType = typename TInputImage::SizeType;
77 using PixelType = typename TInputImage::PixelType;
78 using OutputImagePointer = typename OutputImageType::Pointer;
79 using RegionType = typename OutputImageType::RegionType;
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
100
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
115 itkConceptMacro(OutputHasNumericTraitsCheck,
117
118protected:
120 ~GradientVectorFlowImageFilter() override = default;
121 void
122 PrintSelf(std::ostream & os, Indent indent) const override;
123
124 void
125 GenerateData() override;
126
128 void
130
135 void
137
139 void
141
142private:
143 // parameters;
144 double m_TimeStep{}; // the timestep of each
145 // iteration
146 double m_Steps[Superclass::InputImageDimension]{}; // set to be 1 in all
147 // directions in most cases
148 double m_NoiseLevel{}; // the noise level of the
149 // image
150 int m_IterationNum{}; // the iteration number
151
154
156 InternalImagePointer m_BImage{}; // store the "b" value for every pixel
157
158 typename Superclass::InputImagePointer m_CImage{}; // store the $c_i$ value for
159 // every pixel
160};
161} // end namespace itk
162
163#ifndef ITK_MANUAL_INSTANTIATION
164# include "itkGradientVectorFlowImageFilter.hxx"
165#endif
166
167#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....