PSCF v1.4.0
cuda/FFT.h
1#ifndef PRDC_CUDA_FFT_H
2#define PRDC_CUDA_FFT_H
3
4/*
5* PSCF - Polymer Self-Consistent Field
6*
7* Copyright 2015 - 2025, The Regents of the University of Minnesota
8* Distributed under the terms of the GNU General Public License.
9*/
10
11#include "RField.h"
12#include "CField.h"
13#include "RFieldDft.h"
14
15#include <pscf/math/IntVec.h>
16#include <util/global.h>
17
18#include <cufft.h>
19
20namespace Pscf {
21namespace Prdc {
22namespace Cuda {
23
24 using namespace Util;
25
36 template <int D>
37 class FFT
38 {
39
40 public:
41
45 FFT();
46
50 virtual ~FFT();
51
57 void setup(IntVec<D> const & meshDimensions);
58
59 // Real <-> Complex transforms
60
76 void forwardTransform(RField<D> const & rField, RFieldDft<D>& kField)
77 const;
78
92 void inverseTransformUnsafe(RFieldDft<D>& kField, RField<D>& rField)
93 const;
94
104 void inverseTransformSafe(RFieldDft<D> const & kField,
105 RField<D>& rField) const;
106
107 // Complex <-> Complex transforms
108
128 void forwardTransform(CField<D> const & rField, CField<D>& kField)
129 const;
130
147 void inverseTransform(CField<D> const & kField, CField<D>& rField)
148 const;
149
153 const IntVec<D>& meshDimensions() const;
154
158 bool isSetup() const;
159
160 // Static function
161
173 static
174 void computeKMesh(IntVec<D> const & rMeshDimensions,
175 IntVec<D> & kMeshDimensions,
176 int & kSize );
177
195 static
196 bool hasImplicitInverse(IntVec<D> const & wavevector,
197 IntVec<D> const & meshDimensions);
198
199 private:
200
202 IntVec<D> meshDimensions_;
203
205 mutable RFieldDft<D> kFieldCopy_;
206
208 int rSize_;
209
211 int kSize_;
212
214 cufftHandle rcfPlan_;
215
217 cufftHandle criPlan_;
218
220 cufftHandle ccPlan_;
221
223 bool isSetup_;
224
228 void makePlans();
229
230 };
231
232 // Declarations of explicit specializations
233
234 template <>
235 void FFT<1>::makePlans();
236
237 template <>
238 void FFT<2>::makePlans();
239
240 template <>
241 void FFT<3>::makePlans();
242
243 // Inline functions
244
245 /*
246 * Return the dimensions of the grid for which this was allocated.
247 */
248 template <int D>
249 inline IntVec<D> const & FFT<D>::meshDimensions() const
250 { return meshDimensions_; }
251
252 /*
253 * Has this FFT object been setup?
254 */
255 template <int D>
256 inline bool FFT<D>::isSetup() const
257 { return isSetup_; }
258
259 /*
260 * Does this wavevector have an implicit inverse?
261 */
262 template <int D>
263 inline
264 bool FFT<D>::hasImplicitInverse(IntVec<D> const & wavevector,
266 {
267 int i = wavevector[D-1];
268 int d = meshDimensions[D-1];
269 UTIL_CHECK(i >= 0 && i < D);
270 if ( (i != 0) && (d - i > d / 2) ) {
271 return true;
272 } else {
273 return false;
274 }
275 }
276
277 // Explicit instantiation declarations
278 extern template class FFT<1>;
279 extern template class FFT<2>;
280 extern template class FFT<3>;
281
282} // Cuda
283} // Prdc
284} // Pscf
285#endif
An IntVec<D, T> is a D-component vector of elements of integer type T.
Definition IntVec.h:27
Field of complex values on a regular mesh, allocated on a GPU device.
Definition cuda/CField.h:31
Fourier transform wrapper for real or complex data.
Definition cuda/FFT.h:38
FFT()
Default constructor.
Definition cuda/FFT.tpp:66
static void computeKMesh(IntVec< D > const &rMeshDimensions, IntVec< D > &kMeshDimensions, int &kSize)
Compute dimensions and size of k-space mesh for DFT of real data.
Definition cuda/FFT.tpp:278
void inverseTransform(CField< D > const &kField, CField< D > &rField) const
Compute inverse (complex-to-complex) discrete Fourier transform.
Definition cuda/FFT.tpp:246
bool isSetup() const
Has this FFT object been setup?
Definition cuda/FFT.h:256
static bool hasImplicitInverse(IntVec< D > const &wavevector, IntVec< D > const &meshDimensions)
Does a wavevector have an implicit inverse in DFT mesh for real data?
Definition cuda/FFT.h:264
void forwardTransform(RField< D > const &rField, RFieldDft< D > &kField) const
Compute forward (real-to-complex) discrete Fourier transform.
Definition cuda/FFT.tpp:134
void setup(IntVec< D > const &meshDimensions)
Setup grid dimensions, plans and work space.
Definition cuda/FFT.tpp:97
void inverseTransformUnsafe(RFieldDft< D > &kField, RField< D > &rField) const
Compute inverse (complex-to-real) DFT, overwriting the input.
Definition cuda/FFT.tpp:165
const IntVec< D > & meshDimensions() const
Return the dimensions of the grid for which this was allocated.
Definition cuda/FFT.h:249
void inverseTransformSafe(RFieldDft< D > const &kField, RField< D > &rField) const
Compute inverse (complex-to-real) DFT without overwriting input.
Definition cuda/FFT.tpp:191
virtual ~FFT()
Destructor.
Definition cuda/FFT.tpp:80
Discrete Fourier Transform (DFT) of a real field, allocated on a GPU.
Field of real values on a regular mesh, allocated on a GPU device.
Definition cuda/RField.h:33
File containing preprocessor macros for error handling.
#define UTIL_CHECK(condition)
Assertion macro suitable for serial or parallel production code.
Definition global.h:68
Fields, FFTs, and utilities for periodic boundary conditions (CUDA).
Definition CField.cu:12
Periodic fields and crystallography.
Definition complex.cpp:11
PSCF package top-level namespace.