ITK  6.0.0
Insight Toolkit
itkVnlInverseFFTImageFilter.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 itkVnlInverseFFTImageFilter_h
19#define itkVnlInverseFFTImageFilter_h
20
22
23#include "itkImage.h"
24#include "vnl/algo/vnl_fft_base.h"
25
27
28namespace itk
29{
46template <typename TInputImage,
47 typename TOutputImage = Image<typename TInputImage::PixelType::value_type, TInputImage::ImageDimension>>
48class ITK_TEMPLATE_EXPORT VnlInverseFFTImageFilter : public InverseFFTImageFilter<TInputImage, TOutputImage>
49{
50public:
51 ITK_DISALLOW_COPY_AND_MOVE(VnlInverseFFTImageFilter);
52
54 using InputImageType = TInputImage;
55 using InputPixelType = typename InputImageType::PixelType;
58 using OutputImageType = TOutputImage;
59 using OutputPixelType = typename OutputImageType::PixelType;
61
66
68 itkNewMacro(Self);
69
71 itkOverrideGetNameOfClassMacro(VnlInverseFFTImageFilter);
72
75 static constexpr unsigned int ImageDimension = TOutputImage::ImageDimension;
76 static constexpr unsigned int InputImageDimension = TInputImage::ImageDimension;
77 static constexpr unsigned int OutputImageDimension = TOutputImage::ImageDimension;
78
81
82#ifdef ITK_USE_CONCEPT_CHECKING
83 // Begin concept checking
84 itkConceptMacro(PixelUnsignedIntDivisionOperatorsCheck, (Concept::DivisionOperators<OutputPixelType, unsigned int>));
86 // End concept checking
87#endif
88
89protected:
91 ~VnlInverseFFTImageFilter() override = default;
92
93 void
94 GenerateData() override; // generates output from input
95
96private:
97 using SignalVectorType = vnl_vector<InputPixelType>;
98};
99
100
101// Describe whether input/output are real- or complex-valued
102// for factory registration
103template <>
105{
106 template <typename TUnderlying>
107 using InputPixelType = std::complex<TUnderlying>;
108 template <typename TUnderlying>
109 using OutputPixelType = TUnderlying;
110 using FilterDimensions = std::integer_sequence<unsigned int, 4, 3, 2, 1>;
111};
112
113} // namespace itk
114
115#ifndef ITK_MANUAL_INSTANTIATION
116# include "itkVnlInverseFFTImageFilter.hxx"
117#endif
118
119#endif
Base class for inverse Fast Fourier Transform.
Light weight base class for most itk classes.
VNL-based reverse Fast Fourier Transform.
SizeValueType GetSizeGreatestPrimeFactor() const override
typename InputImageType::SizeValueType InputSizeValueType
~VnlInverseFFTImageFilter() override=default
vnl_vector< InputPixelType > SignalVectorType
typename InputImageType::PixelType InputPixelType
typename OutputImageType::PixelType OutputPixelType
typename InputImageType::SizeType InputSizeType
typename OutputImageType::SizeType OutputSizeType
#define itkConceptMacro(name, concept)
The "itk" namespace contains all Insight Segmentation and Registration Toolkit (ITK) classes....
unsigned long SizeValueType
Definition: itkIntTypes.h:86
std::integer_sequence< unsigned int, 4, 3, 2, 1 > FilterDimensions
Helper defining pixel traits for templated FFT image filters.