forked from InsightSoftwareConsortium/ITKVkFFTBackend
/
itkVkComplexToComplex1DFFTImageFilter.hxx
131 lines (115 loc) · 5.13 KB
/
itkVkComplexToComplex1DFFTImageFilter.hxx
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
/*=========================================================================
*
* 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
*
* http://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 itkVkComplexToComplex1DFFTImageFilter_hxx
#define itkVkComplexToComplex1DFFTImageFilter_hxx
#include "itkVkComplexToComplex1DFFTImageFilter.h"
#include "itkVkGlobalConfiguration.h"
#include "vkFFT.h"
#include "itkImageRegionIterator.h"
#include "itkIndent.h"
#include "itkMetaDataObject.h"
#include "itkProgressReporter.h"
namespace itk
{
template <typename TInputImage, typename TOutputImage>
void
VkComplexToComplex1DFFTImageFilter<TInputImage, TOutputImage>::GenerateData()
{
// get pointers to the input and output
const InputImageType * const input{ this->GetInput() };
OutputImageType * const output{ this->GetOutput() };
if (!input || !output)
{
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
const ProgressReporter progress(this, 0, 1);
// allocate output buffer memory
output->SetBufferedRegion(output->GetRequestedRegion());
output->Allocate();
const SizeType & inputSize{ input->GetLargestPossibleRegion().GetSize() };
const InputPixelType * const inputCPUBuffer{ input->GetBufferPointer() };
OutputPixelType * const outputCPUBuffer{ output->GetBufferPointer() };
itkAssertOrThrowMacro(inputCPUBuffer != nullptr, "No CPU input buffer");
itkAssertOrThrowMacro(outputCPUBuffer != nullptr, "No CPU output buffer");
const SizeValueType inBytes{ input->GetLargestPossibleRegion().GetNumberOfPixels() * sizeof(InputPixelType) };
const SizeValueType outBytes{ output->GetLargestPossibleRegion().GetNumberOfPixels() * sizeof(OutputPixelType) };
itkAssertOrThrowMacro(inBytes == outBytes, "CPU input and output buffers are of different sizes.");
// Mostly use defaults for VkCommon::VkGPU
typename VkCommon::VkGPU vkGPU;
vkGPU.device_id = this->GetDeviceID();
// Describe this filter in VkCommon::VkParameters
typename VkCommon::VkParameters vkParameters;
if (ImageDimension > 0)
vkParameters.X = inputSize[0];
if (ImageDimension > 1)
vkParameters.Y = inputSize[1];
if (ImageDimension > 2)
vkParameters.Z = inputSize[2];
if (std::is_same<RealType, float>::value)
vkParameters.P = VkCommon::PrecisionEnum::FLOAT;
else if (std::is_same<RealType, double>::value)
vkParameters.P = VkCommon::PrecisionEnum::DOUBLE;
else
itkAssertOrThrowMacro(false, "Unsupported type for real numbers.");
vkParameters.fft = VkCommon::FFTEnum::C2C;
vkParameters.PSize = sizeof(RealType);
vkParameters.I = this->GetTransformDirection() == Superclass::TransformDirectionType::INVERSE
? VkCommon::DirectionEnum::INVERSE
: VkCommon::DirectionEnum::FORWARD;
vkParameters.normalized = vkParameters.I == VkCommon::DirectionEnum::INVERSE
? VkCommon::NormalizationEnum::NORMALIZED
: VkCommon::NormalizationEnum::UNNORMALIZED;
for (size_t dim{ 0 }; dim < ImageDimension; ++dim)
{
if (this->GetDirection() != dim)
{
vkParameters.omitDimension[dim] = 1; // omit dimensions other than in the given direction.
}
}
vkParameters.inputCPUBuffer = inputCPUBuffer;
vkParameters.inputBufferBytes = inBytes;
vkParameters.outputCPUBuffer = outputCPUBuffer;
vkParameters.outputBufferBytes = outBytes;
const VkFFTResult resFFT{ m_VkCommon.Run(vkGPU, vkParameters) };
if (resFFT != VKFFT_SUCCESS)
{
std::ostringstream mesg;
mesg << "VkFFT third-party library failed with error code " << resFFT << ".";
itkAssertOrThrowMacro(false, mesg.str());
}
}
template <typename TInputImage, typename TOutputImage>
void
VkComplexToComplex1DFFTImageFilter<TInputImage, TOutputImage>::PrintSelf(std::ostream & os, Indent indent) const
{
Superclass::PrintSelf(os, indent);
os << indent << "UseVkGlobalConfiguration: " << m_UseVkGlobalConfiguration << std::endl;
os << indent << "Local DeviceID: " << m_DeviceID << std::endl;
os << indent << "Global DeviceID: " << VkGlobalConfiguration::GetDeviceID() << std::endl;
os << indent << "Preferred DeviceID: " << this->GetDeviceID() << std::endl;
}
template <typename TInputImage, typename TOutputImage>
typename VkComplexToComplex1DFFTImageFilter<TInputImage, TOutputImage>::SizeValueType
VkComplexToComplex1DFFTImageFilter<TInputImage, TOutputImage>::GetSizeGreatestPrimeFactor() const
{
return SizeValueType{ m_VkCommon.GetGreatestPrimeFactor() };
}
} // end namespace itk
#endif // _itkVkComplexToComplex1DFFTImageFilter_hxx