File: itkFFTWComplexToComplex1DFFTImageFilter.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 (204 lines) | stat: -rw-r--r-- 6,838 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
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
/*=========================================================================
 *
 *  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 itkFFTWComplexToComplex1DFFTImageFilter_hxx
#define itkFFTWComplexToComplex1DFFTImageFilter_hxx

#include "itkComplexToComplex1DFFTImageFilter.hxx"

#include "itkFFTWCommonExtended.h"
#include "itkIndent.h"
#include "itkImageLinearConstIteratorWithIndex.h"
#include "itkImageLinearIteratorWithIndex.h"
#include "itkMetaDataObject.h"

#if defined(ITK_USE_FFTWF) || defined(ITK_USE_FFTWD)

namespace itk
{

template <typename TInputImage, typename TOutputImage>
FFTWComplexToComplex1DFFTImageFilter<TInputImage, TOutputImage>::FFTWComplexToComplex1DFFTImageFilter()

{
  // We cannot split over the FFT direction
  this->m_ImageRegionSplitter = ImageRegionSplitterDirection::New();
  this->DynamicMultiThreadingOff();
}


template <typename TInputImage, typename TOutputImage>
FFTWComplexToComplex1DFFTImageFilter<TInputImage, TOutputImage>::~FFTWComplexToComplex1DFFTImageFilter()
{
  if (m_PlanComputed)
  {
    this->DestroyPlans();
  }
}


template <typename TInputImage, typename TOutputImage>
void
FFTWComplexToComplex1DFFTImageFilter<TInputImage, TOutputImage>::DestroyPlans()
{
  for (unsigned int i = 0; i < m_PlanArray.size(); i++)
  {
    FFTW1DProxyType::DestroyPlan(m_PlanArray[i]);
    delete[] m_InputBufferArray[i];
    delete[] m_OutputBufferArray[i];
    this->m_PlanComputed = false;
  }
}


template <typename TInputImage, typename TOutputImage>
const ImageRegionSplitterBase *
FFTWComplexToComplex1DFFTImageFilter<TInputImage, TOutputImage>::GetImageRegionSplitter() const
{
  return this->m_ImageRegionSplitter.GetPointer();
}


template <typename TInputImage, typename TOutputImage>
void
FFTWComplexToComplex1DFFTImageFilter<TInputImage, TOutputImage>::BeforeThreadedGenerateData()
{
  Superclass::BeforeThreadedGenerateData();

  this->m_ImageRegionSplitter->SetDirection(this->GetDirection());

  OutputImageType * outputPtr = this->GetOutput();

  const typename OutputImageType::SizeType & outputSize = outputPtr->GetRequestedRegion().GetSize();
  const unsigned int                         lineSize = outputSize[this->m_Direction];

  if (this->m_PlanComputed)
  {
    // if the image sizes aren't the same,
    // we have to compute the plan again
    if (this->m_LastImageSize != lineSize)
    {
      this->DestroyPlans();
    }
  }
  if (!this->m_PlanComputed)
  {
    const int threads = this->GetNumberOfWorkUnits();
    m_PlanArray.resize(threads);
    m_InputBufferArray.resize(threads);
    m_OutputBufferArray.resize(threads);
    for (int i = 0; i < threads; i++)
    {
      try
      {
        m_InputBufferArray[i] = new typename FFTW1DProxyType::ComplexType[lineSize];
        m_OutputBufferArray[i] = new typename FFTW1DProxyType::ComplexType[lineSize];
      }
      catch (const std::bad_alloc &)
      {
        itkExceptionMacro("Problem allocating memory for internal computations");
      }
      if (this->m_TransformDirection == Superclass::DIRECT)
      {
        m_PlanArray[i] = FFTW1DProxyType::Plan_dft_1d(
          lineSize, m_InputBufferArray[i], m_OutputBufferArray[i], FFTW_FORWARD, FFTW_ESTIMATE, 1);
      }
      else // m_TransformDirection == INVERSE
      {
        m_PlanArray[i] = FFTW1DProxyType::Plan_dft_1d(
          lineSize, m_InputBufferArray[i], m_OutputBufferArray[i], FFTW_BACKWARD, FFTW_ESTIMATE, 1);
      }
    }
    this->m_LastImageSize = lineSize;
    this->m_PlanComputed = true;
  }
}


template <typename TInputImage, typename TOutputImage>
void
FFTWComplexToComplex1DFFTImageFilter<TInputImage, TOutputImage>::ThreadedGenerateData(
  const OutputImageRegionType & outputRegion,
  ThreadIdType                  threadID)
{
  // get pointers to the input and output
  const InputImageType * inputPtr = this->GetInput();
  OutputImageType *      outputPtr = this->GetOutput();

  const typename OutputImageType::SizeType & outputSize = outputPtr->GetRequestedRegion().GetSize();
  const unsigned int                         lineSize = outputSize[this->m_Direction];

  using InputIteratorType = itk::ImageLinearConstIteratorWithIndex<InputImageType>;
  using OutputIteratorType = itk::ImageLinearIteratorWithIndex<OutputImageType>;
  InputIteratorType inputIt(inputPtr, outputRegion);
  // the output region should be the same as the input region in the non-fft directions
  OutputIteratorType outputIt(outputPtr, outputRegion);

  inputIt.SetDirection(this->m_Direction);
  outputIt.SetDirection(this->m_Direction);

  typename InputIteratorType::PixelType *  inputBufferIt;
  typename OutputIteratorType::PixelType * outputBufferIt;

  // for every fft line
  for (inputIt.GoToBegin(), outputIt.GoToBegin(); !inputIt.IsAtEnd(); outputIt.NextLine(), inputIt.NextLine())
  {
    // copy the input line into our buffer
    inputIt.GoToBeginOfLine();
    inputBufferIt = reinterpret_cast<typename InputIteratorType::PixelType *>(m_InputBufferArray[threadID]);
    while (!inputIt.IsAtEndOfLine())
    {
      *inputBufferIt = inputIt.Get();
      ++inputIt;
      ++inputBufferIt;
    }

    // do the transform
    FFTW1DProxyType::Execute(m_PlanArray[threadID]);

    if (this->m_TransformDirection == Superclass::DIRECT)
    {
      // copy the output from the buffer into our line
      outputBufferIt = reinterpret_cast<typename OutputIteratorType::PixelType *>(m_OutputBufferArray[threadID]);
      outputIt.GoToBeginOfLine();
      while (!outputIt.IsAtEndOfLine())
      {
        outputIt.Set(*outputBufferIt);
        ++outputIt;
        ++outputBufferIt;
      }
    }
    else // m_TransformDirection == INVERSE
    {
      // copy the output from the buffer into our line
      outputBufferIt = reinterpret_cast<typename OutputIteratorType::PixelType *>(m_OutputBufferArray[threadID]);
      outputIt.GoToBeginOfLine();
      while (!outputIt.IsAtEndOfLine())
      {
        outputIt.Set(*outputBufferIt / static_cast<typename OutputIteratorType::PixelType>(lineSize));
        ++outputIt;
        ++outputBufferIt;
      }
    }
  }
}

} // namespace itk

#endif // defined( ITK_USE_FFTWF ) || defined( ITK_USE_FFTWD )

#endif //_itkFFTWComplexToComplex1DFFTImageFilter_hxx