ITK 6.0.0
Insight Toolkit
 
Loading...
Searching...
No Matches
itkFFTDiscreteGaussianImageFilter.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 itkFFTDiscreteGaussianImageFilter_h
19#define itkFFTDiscreteGaussianImageFilter_h
20
23#include "itkImage.h"
24#include "itkMacro.h"
25#include "ITKSmoothingExport.h"
26
27namespace itk
28{
34{
35public:
51 enum class KernelSource : uint8_t
52 {
55 };
56};
57// Define how to print enumeration
58extern ITKSmoothing_EXPORT std::ostream &
60
80
81template <typename TInputImage, typename TOutputImage = TInputImage>
82class ITK_TEMPLATE_EXPORT FFTDiscreteGaussianImageFilter : public DiscreteGaussianImageFilter<TInputImage, TOutputImage>
83{
84public:
85 ITK_DISALLOW_COPY_AND_MOVE(FFTDiscreteGaussianImageFilter);
86
92
94 itkNewMacro(Self);
95
97 itkOverrideGetNameOfClassMacro(FFTDiscreteGaussianImageFilter);
98
102
109
113
116 static constexpr unsigned int ImageDimension = Superclass::ImageDimension;
117
119 using typename Superclass::RealOutputPixelType;
120 using typename Superclass::RealOutputImageType;
124
130
132 using typename Superclass::ArrayType;
133 using typename Superclass::SigmaArrayType;
134 using typename Superclass::ScalarRealType;
135
137 using typename Superclass::KernelType;
138 using typename Superclass::RadiusType;
139
142
144
145 void
147
150
151 itkGetConstMacro(KernelImage, typename RealImageType::Pointer);
152
153protected:
156
158 void
160
166 void
167 GenerateData() override;
168
171 auto
173
174private:
175 /* Enum for choosing among available kernel generation methods */
178
179 /* Kernel image is allocated with GenerateKernelImage() */
181
182 /* Persist mini-pipeline filter to minimize construction costs
183 * on repeated calls to GenerateData() */
185};
186} // end namespace itk
187
188#ifndef ITK_MANUAL_INSTANTIATION
189# include "itkFFTDiscreteGaussianImageFilter.hxx"
190#endif
191
192#endif
ZeroFluxNeumannBoundaryCondition< TInputImage > InputDefaultBoundaryConditionType
typename NumericTraits< InputPixelType >::ValueType InputPixelValueType
GaussianOperator< RealOutputPixelValueType, ImageDimension > KernelType
ImageBoundaryCondition< RealOutputImageType > * RealBoundaryConditionPointerType
ZeroFluxNeumannBoundaryCondition< RealOutputImageType > RealDefaultBoundaryConditionType
typename NumericTraits< RealOutputPixelType >::ValueType RealOutputPixelValueType
typename NumericTraits< OutputPixelType >::ValueType OutputPixelValueType
typename NumericTraits< OutputPixelType >::RealType RealOutputPixelType
Convolve a given image with an arbitrary image kernel using multiplication in the Fourier domain.
Contains all enum classes used by FFTDiscreteGaussianImageFilter class.
FFTConvolutionImageFilter< RealImageType, RealImageType, OutputImageType > ConvolutionImageFilterType
void GenerateInputRequestedRegion() override
typename Superclass::OutputInternalPixelType OutputInternalPixelType
typename Superclass::InputPixelType InputPixelType
FFTDiscreteGaussianImageFilterEnums::KernelSource m_KernelSource
typename Superclass::OutputPixelType OutputPixelType
typename Superclass::OutputPixelValueType OutputPixelValueType
ConvolutionImageFilterType::Pointer m_ConvolutionImageFilter
typename Superclass::OutputImageType OutputImageType
typename Superclass::InputInternalPixelType InputInternalPixelType
typename Superclass::InputImageType InputImageType
typename Superclass::InputPixelValueType InputPixelValueType
~FFTDiscreteGaussianImageFilter() override=default
auto GenerateKernelImage() -> RealImageType *
void SetInputBoundaryCondition(BoundaryConditionType *) override
DiscreteGaussianImageFilter< TInputImage, TOutputImage > Superclass
Implements transparent reference counting.
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
ITKCommon_EXPORT std::ostream & operator<<(std::ostream &out, typename AnatomicalOrientation::CoordinateEnum value)