File: itkInverseFFTImageFilter.h

package info (click to toggle)
insighttoolkit5 5.4.4-1
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 704,404 kB
  • sloc: cpp: 783,697; ansic: 628,724; xml: 44,704; fortran: 34,250; python: 22,874; sh: 4,078; pascal: 2,636; lisp: 2,158; makefile: 461; yacc: 328; asm: 205; perl: 203; lex: 146; tcl: 132; javascript: 98; csh: 81
file content (109 lines) | stat: -rw-r--r-- 3,674 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
/*=========================================================================
 *
 *  Copyright NumFOCUS
 *
 *  Licensed under the Apache License, Version 2.0 (the "License");
 *  you may not use this file except in compliance with the License.
 *  You may obtain a copy of the License at
 *
 *         https://www.apache.org/licenses/LICENSE-2.0.txt
 *
 *  Unless required by applicable law or agreed to in writing, software
 *  distributed under the License is distributed on an "AS IS" BASIS,
 *  WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
 *  See the License for the specific language governing permissions and
 *  limitations under the License.
 *
 *=========================================================================*/
#ifndef itkInverseFFTImageFilter_h
#define itkInverseFFTImageFilter_h

#include "itkImageToImageFilter.h"
#include "itkMacro.h"

namespace itk
{
/**
 * \class InverseFFTImageFilter
 *
 * \brief Base class for inverse Fast Fourier Transform.
 *
 * This is a base class for the "inverse" or "reverse" Discrete Fourier
 * Transform.  This is an abstract base class: the actual implementation is
 * provided by the best child available on the system when the object is
 * created via the object factory system.
 *
 * This class transforms a full complex image with Hermitian symmetry into
 * its real spatial domain representation.  If the input does not have
 * Hermitian symmetry, the imaginary component is discarded.
 *
 * \ingroup FourierTransform
 *
 * \sa ForwardFFTImageFilter, InverseFFTImageFilter
 * \ingroup ITKFFT
 *
 * \sphinx
 * \sphinxexample{Filtering/FFT/ComputeInverseFFTOfImage,Compute Inverse FFT Of Image}
 * \endsphinx
 */
template <typename TInputImage,
          typename TOutputImage = Image<typename TInputImage::PixelType::value_type, TInputImage::ImageDimension>>
class ITK_TEMPLATE_EXPORT InverseFFTImageFilter : public ImageToImageFilter<TInputImage, TOutputImage>

{
public:
  ITK_DISALLOW_COPY_AND_MOVE(InverseFFTImageFilter);

  /** Standard class type aliases. */
  using InputImageType = TInputImage;
  using InputPixelType = typename InputImageType::PixelType;
  using OutputImageType = TOutputImage;
  using OutputPixelType = typename OutputImageType::PixelType;

  using Self = InverseFFTImageFilter;
  using Superclass = ImageToImageFilter<InputImageType, OutputImageType>;
  using Pointer = SmartPointer<Self>;
  using ConstPointer = SmartPointer<const Self>;

  static constexpr unsigned int ImageDimension = InputImageType::ImageDimension;

  /** \see LightObject::GetNameOfClass() */
  itkOverrideGetNameOfClassMacro(InverseFFTImageFilter);

  /** Customized object creation methods that support configuration-based
   * selection of FFT implementation.
   *
   * Default implementation is VnlFFT. */
  itkFactoryOnlyNewMacro(Self);

  /* Return the preferred greatest prime factor supported for the input image
   * size. Defaults to 2 as many implementations work only for sizes that are
   * power of 2.
   */
  virtual SizeValueType
  GetSizeGreatestPrimeFactor() const;

protected:
  InverseFFTImageFilter() = default;
  ~InverseFFTImageFilter() override = default;

  /** This class requires the entire input. */
  void
  GenerateInputRequestedRegion() override;

  /** Sets the output requested region to the largest possible output
   * region. */
  void
  EnlargeOutputRequestedRegion(DataObject * itkNotUsed(output)) override;
};
} // end namespace itk

#ifndef ITK_MANUAL_INSTANTIATION
#  include "itkInverseFFTImageFilter.hxx"
#endif

#ifdef ITK_FFTIMAGEFILTERINIT_FACTORY_REGISTER_MANAGER
#  include "itkFFTImageFilterInitFactoryRegisterManager.h"
#endif

#endif