File: itkFFTWInverseFFTImageFilter.hxx

package info (click to toggle)
insighttoolkit5 5.4.3-5
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid, trixie
  • size: 704,384 kB
  • sloc: cpp: 783,592; ansic: 628,724; xml: 44,704; fortran: 34,250; python: 22,874; sh: 4,078; pascal: 2,636; lisp: 2,158; makefile: 464; yacc: 328; asm: 205; perl: 203; lex: 146; tcl: 132; javascript: 98; csh: 81
file content (137 lines) | stat: -rw-r--r-- 4,661 bytes parent folder | download
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
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
/*=========================================================================
 *
 *  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 itkFFTWInverseFFTImageFilter_hxx
#define itkFFTWInverseFFTImageFilter_hxx

#include "itkFullToHalfHermitianImageFilter.h"

#include "itkImageRegionIterator.h"
#include "itkProgressReporter.h"
#include "itkMultiThreaderBase.h"

namespace itk
{

template <typename TInputImage, typename TOutputImage>
FFTWInverseFFTImageFilter<TInputImage, TOutputImage>::FFTWInverseFFTImageFilter()
{
#ifndef ITK_USE_CUFFTW
  m_PlanRigor = FFTWGlobalConfiguration::GetPlanRigor();
#endif
  this->DynamicMultiThreadingOn();
}

template <typename TInputImage, typename TOutputImage>
void
FFTWInverseFFTImageFilter<TInputImage, TOutputImage>::BeforeThreadedGenerateData()
{
  // Get pointers to the input and output.
  typename InputImageType::ConstPointer inputPtr = this->GetInput();
  typename OutputImageType::Pointer     outputPtr = this->GetOutput();

  if (!inputPtr || !outputPtr)
  {
    return;
  }

  // We don't have a nice progress to report, but at least this simple line
  // reports the beginning and the end of the process.
  ProgressReporter progress(this, 0, 1);

  // Allocate output buffer memory.
  outputPtr->SetBufferedRegion(outputPtr->GetRequestedRegion());
  outputPtr->Allocate();

  const InputSizeType  inputSize = inputPtr->GetLargestPossibleRegion().GetSize();
  const OutputSizeType outputSize = outputPtr->GetLargestPossibleRegion().GetSize();

  // Figure out sizes.
  // Size of input and output aren't the same which is handled in the superclass,
  // sort of.
  // The input size and output size only differ in the fastest moving dimension.
  unsigned int totalOutputSize = 1;
  unsigned int totalInputSize = 1;

  for (unsigned int i = 0; i < ImageDimension; ++i)
  {
    totalOutputSize *= outputSize[i];
    totalInputSize *= inputSize[i];
  }

  // Cut the full complex image to just the portion needed by FFTW.
  using FullToHalfFilterType = FullToHalfHermitianImageFilter<InputImageType>;
  auto fullToHalfFilter = FullToHalfFilterType::New();
  fullToHalfFilter->SetInput(this->GetInput());
  fullToHalfFilter->SetNumberOfWorkUnits(this->GetNumberOfWorkUnits());
  fullToHalfFilter->UpdateLargestPossibleRegion();

  auto * in = (typename FFTWProxyType::ComplexType *)fullToHalfFilter->GetOutput()->GetBufferPointer();

  OutputPixelType *                out = outputPtr->GetBufferPointer();
  typename FFTWProxyType::PlanType plan;

  int sizes[ImageDimension];
  for (unsigned int i = 0; i < ImageDimension; ++i)
  {
    sizes[(ImageDimension - 1) - i] = outputSize[i];
  }

  plan = FFTWProxyType::Plan_dft_c2r(
    ImageDimension, sizes, in, out, m_PlanRigor, MultiThreaderBase::GetGlobalDefaultNumberOfThreads(), false);
  FFTWProxyType::Execute(plan);

  // Some cleanup.
  FFTWProxyType::DestroyPlan(plan);
}

template <typename TInputImage, typename TOutputImage>
void
FFTWInverseFFTImageFilter<TInputImage, TOutputImage>::DynamicThreadedGenerateData(
  const OutputImageRegionType & outputRegionForThread)
{
  using IteratorType = ImageRegionIterator<OutputImageType>;
  unsigned long totalOutputSize = this->GetOutput()->GetRequestedRegion().GetNumberOfPixels();
  IteratorType  it(this->GetOutput(), outputRegionForThread);
  while (!it.IsAtEnd())
  {
    it.Set(it.Value() / totalOutputSize);
    ++it;
  }
}

template <typename TInputImage, typename TOutputImage>
void
FFTWInverseFFTImageFilter<TInputImage, TOutputImage>::PrintSelf(std::ostream & os, Indent indent) const
{
  Superclass::PrintSelf(os, indent);

#ifndef ITK_USE_CUFFTW
  os << indent << "PlanRigor: " << FFTWGlobalConfiguration::GetPlanRigorName(m_PlanRigor) << " (" << m_PlanRigor << ')'
     << std::endl;
#endif
}

template <typename TInputImage, typename TOutputImage>
SizeValueType
FFTWInverseFFTImageFilter<TInputImage, TOutputImage>::GetSizeGreatestPrimeFactor() const
{
  return FFTWProxyType::GREATEST_PRIME_FACTOR;
}

} // namespace itk
#endif // _itkFFTWInverseFFTImageFilter_hxx