PSCF v1.3
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#include "types.h"
15
16#include <pscf/math/IntVec.h>
17#include <util/global.h>
18
19#include <cufft.h>
20
21namespace Pscf {
22namespace Prdc {
23namespace Cuda {
24
25 using namespace Util;
26
37 template <int D>
38 class FFT
39 {
40
41 public:
42
46 FFT();
47
51 virtual ~FFT();
52
59
60 // Real <-> Complex transforms
61
77 void forwardTransform(RField<D> const & rField, RFieldDft<D>& kField)
78 const;
79
94 const;
95
105 void inverseTransformSafe(RFieldDft<D> const & kField,
106 RField<D>& rField) const;
107
108 // Complex <-> Complex transforms
109
129 void forwardTransform(CField<D> const & rField, CField<D>& kField)
130 const;
131
148 void inverseTransform(CField<D> const & kField, CField<D>& rField)
149 const;
150
154 const IntVec<D>& meshDimensions() const;
155
159 bool isSetup() const;
160
161 // Static function
162
174 static
175 void computeKMesh(IntVec<D> const & rMeshDimensions,
176 IntVec<D> & kMeshDimensions,
177 int & kSize );
178
196 static
197 bool hasImplicitInverse(IntVec<D> const & wavevector,
198 IntVec<D> const & meshDimensions);
199
200 private:
201
203 IntVec<D> meshDimensions_;
204
206 mutable RFieldDft<D> kFieldCopy_;
207
209 int rSize_;
210
212 int kSize_;
213
215 cufftHandle rcfPlan_;
216
218 cufftHandle criPlan_;
219
221 cufftHandle ccPlan_;
222
224 bool isSetup_;
225
229 void makePlans();
230
231 };
232
233 // Declarations of explicit specializations
234
235 template <>
236 void FFT<1>::makePlans();
237
238 template <>
239 void FFT<2>::makePlans();
240
241 template <>
242 void FFT<3>::makePlans();
243
244 // Inline functions
245
246 /*
247 * Return the dimensions of the grid for which this was allocated.
248 */
249 template <int D>
250 inline IntVec<D> const & FFT<D>::meshDimensions() const
251 { return meshDimensions_; }
252
253 /*
254 * Has this FFT object been setup?
255 */
256 template <int D>
257 inline bool FFT<D>::isSetup() const
258 { return isSetup_; }
259
260 /*
261 * Does this wavevector have an implicit inverse?
262 */
263 template <int D>
264 inline
265 bool FFT<D>::hasImplicitInverse(IntVec<D> const & wavevector,
266 IntVec<D> const & meshDimensions)
267 {
268 int i = wavevector[D-1];
269 int d = meshDimensions[D-1];
270 if ((i != 0) && (d - i > d/2 + 1)) {
271 return true;
272 } else {
273 return false;
274 }
275 }
276
277 #ifndef PRDC_CUDA_FFT_TPP
278 // Suppress implicit instantiation
279 extern template class FFT<1>;
280 extern template class FFT<2>;
281 extern template class FFT<3>;
282 #endif
283
284} // Cuda
285} // Prdc
286} // Pscf
287#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
Fourier transform wrapper.
Definition cpu/FFT.h:38
Fourier transform of a real field on an FFT mesh.
Field of real double precision values on an FFT mesh.
Definition cpu/RField.h:29
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.
FFT()
Default constructor.
void inverseTransform(CField< D > const &kField, CField< D > &rField) const
Compute inverse (complex-to-complex) discrete Fourier transform.
virtual ~FFT()
Destructor.
bool isSetup() const
Has this FFT object been setup?
const IntVec< D > & meshDimensions() const
Return the dimensions of the grid for which this was allocated.
void forwardTransform(RField< D > const &rField, RFieldDft< D > &kField) const
Compute forward (real-to-complex) discrete Fourier transform.
void setup(IntVec< D > const &meshDimensions)
Setup grid dimensions, plans and work space.
void inverseTransformUnsafe(RFieldDft< D > &kField, RField< D > &rField) const
Compute inverse (complex-to-real) DFT, overwriting the input.
void inverseTransformSafe(RFieldDft< D > const &kField, RField< D > &rField) const
Compute inverse (complex-to-real) DFT without overwriting input.
void forwardTransform(CField< D > const &rField, CField< D > &kField) const
Compute forward (complex-to-complex) discrete Fourier transform.
static bool hasImplicitInverse(IntVec< D > const &wavevector, IntVec< D > const &meshDimensions)
Does a wavevector have an implicit inverse in DFT for real data?
File containing preprocessor macros for error handling.
Fields, FFTs, and utilities for periodic boundary conditions (CUDA)
Definition Reduce.cpp:14
Periodic fields and crystallography.
Definition CField.cpp:11
PSCF package top-level namespace.
Definition param_pc.dox:1