PSCF v1.2
cpu/FFT.tpp
1#ifndef PRDC_CPU_FFT_TPP
2#define PRDC_CPU_FFT_TPP
3
4/*
5* PSCF Package
6*
7* Copyright 2016 - 2022, The Regents of the University of Minnesota
8* Distributed under the terms of the GNU General Public License.
9*/
10
11#include "FFT.h"
12
13/*
14* A note about const_casts:
15*
16* The FFTW library is used in this file to perform discrete Fourier
17* transforms. FFTW's complex-to-real inverse transform overwrites its
18* input array, but all other out-of-place transforms leave the input
19* array unaltered. However, all transforms in the FFTW library require
20* non-const pointers to the input array, even though they do not alter
21* the array.
22*
23* In order to maintain const-correctness in PSCF, this FFT class accepts
24* const input arrays for its methods that perform a Fourier transform,
25* unless the transform is expected to modify / overwrite its input (as
26* is the case for complex-to-real inverse transforms). This denotes to
27* the caller of the method that the input array will not be altered,
28* which is an accurate representation of the expected behavior.
29*
30* However, the const-correctness of this FFT class creates a conflict
31* with the FFTW library's requirement of non-const inputs. This conflict
32* is resolved using a const_cast, in which the const pointer to the input
33* array is made non-const when passed into FFTW functions. The use of
34* const_cast is reserved only for these few select cases in which we are
35* confident that the input array will not be modified.
36*
37* For more information about the relevant FFTW methods, see the FFTW
38* documentation at https://www.fftw.org/fftw3_doc/index.html. In Section
39* 4.3.2, it is stated that, by default, "an out-of-place transform must
40* not change its input array," except for complex-to-real transforms, in
41* which case "no input-preserving algorithms are implemented." Finally,
42* we note that the unit tests for this FFT class check that the input
43* array is unaltered, allowing developers to continually ensure that
44* the FFTW functions do not modify their input unexpectedly.
45*/
46
47namespace Pscf {
48namespace Prdc {
49namespace Cpu {
50
51 using namespace Util;
52
53 /*
54 * Default constructor.
55 */
56 template <int D>
58 : meshDimensions_(0),
59 rSize_(0),
60 kSize_(0),
61 rcfPlan_(0),
62 criPlan_(0),
63 ccfPlan_(0),
64 cciPlan_(0),
65 isSetup_(false)
66 {}
67
68 /*
69 * Destructor.
70 */
71 template <int D>
73 {
74 if (rcfPlan_) {
75 fftw_destroy_plan(rcfPlan_);
76 }
77 if (criPlan_) {
78 fftw_destroy_plan(criPlan_);
79 }
80 if (ccfPlan_) {
81 fftw_destroy_plan(ccfPlan_);
82 }
83 if (cciPlan_) {
84 fftw_destroy_plan(cciPlan_);
85 }
86 }
87
88 /*
89 * Setup mesh dimensions, work memory and FFT plans.
90 */
91 template <int D>
92 void FFT<D>::setup(IntVec<D> const & meshDimensions)
93 {
94 // Precondition
95 UTIL_CHECK(!isSetup_);
96
97 // Set and check mesh dimensions
98 rSize_ = 1;
99 kSize_ = 1;
100 for (int i = 0; i < D; ++i) {
101 UTIL_CHECK(meshDimensions[i] > 0);
102 meshDimensions_[i] = meshDimensions[i];
103 rSize_ *= meshDimensions[i];
104 if (i < D - 1) {
105 kSize_ *= meshDimensions[i];
106 } else {
107 kSize_ *= (meshDimensions[i]/2 + 1);
108 }
109 }
110
111 // Allocate kFieldCopy_ array if necessary
112 if (!kFieldCopy_.isAllocated()) {
113 kFieldCopy_.allocate(meshDimensions);
114 } else {
115 if (kFieldCopy_.capacity() != kSize_) {
116 kFieldCopy_.deallocate();
117 kFieldCopy_.allocate(meshDimensions);
118 }
119 }
120 UTIL_CHECK(meshDimensions == kFieldCopy_.meshDimensions());
121 UTIL_CHECK(kFieldCopy_.capacity() == kSize_);
122
123 // Create temporary RField and CField objects used for plans
124 RField<D> rField;
125 rField.allocate(meshDimensions);
126 UTIL_CHECK(meshDimensions == rField.meshDimensions());
127 UTIL_CHECK(rField.capacity() == rSize_);
128
129 CField<D> cFieldIn;
130 cFieldIn.allocate(meshDimensions);
131 UTIL_CHECK(meshDimensions == cFieldIn.meshDimensions());
132 UTIL_CHECK(cFieldIn.capacity() == rSize_);
133
134 CField<D> cFieldOut;
135 cFieldOut.allocate(meshDimensions);
136 UTIL_CHECK(meshDimensions == cFieldOut.meshDimensions());
137 UTIL_CHECK(cFieldOut.capacity() == rSize_);
138
139 // Make FFTW plans (see explicit specializations in FFT.cpp)
140 makePlans(rField, kFieldCopy_, cFieldIn, cFieldOut);
141
142 isSetup_ = true;
143 }
144
145 // Real <-> Complex Transforms
146
147 /*
148 * Execute real-to-complex forward transform.
149 */
150 template <int D>
152 RFieldDft<D>& kField)
153 const
154 {
155 UTIL_CHECK(isSetup_)
156 UTIL_CHECK(rField.capacity() == rSize_);
157 UTIL_CHECK(kField.capacity() == kSize_);
158 UTIL_CHECK(rField.meshDimensions() == meshDimensions_);
159 UTIL_CHECK(kField.meshDimensions() == meshDimensions_);
160
161 // Execute preplanned forward transform
162 // (See note at top of file explaining this use of const_cast)
163 fftw_execute_dft_r2c(rcfPlan_, const_cast<double*>(&rField[0]),
164 &kField[0]);
165
166 // Rescale the resulting array
167 double scale = 1.0/double(rSize_);
168 for (int i = 0; i < kSize_; ++i) {
169 kField[i][0] *= scale;
170 kField[i][1] *= scale;
171 }
172 }
173
174 /*
175 * Execute inverse (complex-to-real) transform.
176 */
177 template <int D>
178 void
180 const
181 {
182 UTIL_CHECK(isSetup_)
183 UTIL_CHECK(rField.capacity() == rSize_);
184 UTIL_CHECK(kField.capacity() == kSize_);
185 UTIL_CHECK(rField.meshDimensions() == meshDimensions_);
186 UTIL_CHECK(kField.meshDimensions() == meshDimensions_);
187
188 // Execute preplanned inverse transform
189 fftw_execute_dft_c2r(criPlan_, &kField[0], &rField[0]);
190 }
191
192 /*
193 * Execute inverse (complex-to-real) transform without destroying input.
194 */
195 template <int D>
196 void
198 RField<D>& rField)
199 const
200 {
201 UTIL_CHECK(kFieldCopy_.capacity() == kField.capacity());
202
203 kFieldCopy_ = kField;
204 inverseTransformUnsafe(kFieldCopy_, rField);
205 }
206
207 // Complex <-> Complex Transforms
208
209 /*
210 * Execute complex-to-complex forward transform.
211 */
212 template <int D>
213 void
215 const
216 {
217 UTIL_CHECK(isSetup_)
218 UTIL_CHECK(rField.capacity() == rSize_);
219 UTIL_CHECK(rField.meshDimensions() == meshDimensions_);
220 UTIL_CHECK(kField.capacity() == rSize_);
221 UTIL_CHECK(kField.meshDimensions() == meshDimensions_);
222
223 // Execute preplanned forward transform
224 // (See note at top of file explaining this use of const_cast)
225 fftw_execute_dft(ccfPlan_, const_cast<fftw_complex*>(&rField[0]),
226 &kField[0]);
227
228 // Rescale the resulting array
229 double scale = 1.0/double(rSize_);
230 for (int i = 0; i < rSize_; ++i) {
231 kField[i][0] *= scale;
232 kField[i][1] *= scale;
233 }
234 }
235
236 /*
237 * Execute inverse (complex-to-complex) transform.
238 */
239 template <int D>
240 void
242 const
243 {
244 UTIL_CHECK(isSetup_)
245 UTIL_CHECK(rField.capacity() == rSize_);
246 UTIL_CHECK(kField.capacity() == rSize_);
247 UTIL_CHECK(rField.meshDimensions() == meshDimensions_);
248 UTIL_CHECK(kField.meshDimensions() == meshDimensions_);
249
250 // Execute preplanned inverse transform
251 // (See note at top of file explaining this use of const_cast)
252 fftw_execute_dft(cciPlan_, const_cast<fftw_complex*>(&kField[0]),
253 &rField[0]);
254 }
255
256}
257}
258}
259#endif
An IntVec<D, T> is a D-component vector of elements of integer type T.
Definition IntVec.h:27
Field of complex double precision values on an FFT mesh.
Definition cpu/CField.h:30
const IntVec< D > & meshDimensions() const
Return mesh dimensions by constant reference.
Definition cpu/CField.h:127
void allocate(const IntVec< D > &meshDimensions)
Allocate the underlying C array for an FFT grid.
void inverseTransformUnsafe(RFieldDft< D > &in, RField< D > &out) const
Compute inverse (complex-to-real) DFT, overwriting the input.
Definition cpu/FFT.tpp:179
FFT()
Default constructor.
Definition cpu/FFT.tpp:57
void forwardTransform(RField< D > const &in, RFieldDft< D > &out) const
Compute forward (real-to-complex) discrete Fourier transform.
Definition cpu/FFT.tpp:151
void inverseTransform(CField< D > const &in, CField< D > &out) const
Compute complex-to-complex inverse Fourier transform.
Definition cpu/FFT.tpp:241
void setup(IntVec< D > const &meshDimensions)
Setup grid dimensions, FFT plans and work space.
Definition cpu/FFT.tpp:92
virtual ~FFT()
Destructor.
Definition cpu/FFT.tpp:72
void inverseTransformSafe(RFieldDft< D > const &in, RField< D > &out) const
Compute inverse (complex-to-real) DFT without overwriting input.
Definition cpu/FFT.tpp:197
Fourier transform of a real field on an FFT mesh.
IntVec< D > const & meshDimensions() const
Return vector of spatial mesh dimensions by constant reference.
Field of real double precision values on an FFT mesh.
void allocate(IntVec< D > const &meshDimensions)
Allocate the underlying C array for an FFT grid.
const IntVec< D > & meshDimensions() const
Return mesh dimensions by constant reference.
Definition cpu/RField.h:112
int capacity() const
Return allocated size.
Definition Array.h:159
#define UTIL_CHECK(condition)
Assertion macro suitable for serial or parallel production code.
Definition global.h:68
PSCF package top-level namespace.
Definition param_pc.dox:1
Utility classes for scientific computation.