PSCF v1.1
pspg/field/FFT.h
1#ifndef PSPG_FFT_H
2#define PSPG_FFT_H
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 <pspg/field/RDField.h>
12#include <pspg/field/RDFieldDft.h>
13#include <pspg/math/GpuResources.h>
14#include <pscf/math/IntVec.h>
15#include <util/global.h>
16
17#include <cufft.h>
18#include <cuda.h>
19#include <cuda_runtime.h>
20
21namespace Pscf {
22namespace Pspg {
23
24 using namespace Util;
25 using namespace Pscf;
26
32 template <int D>
33 class FFT
34 {
35
36 public:
37
41 FFT();
42
46 virtual ~FFT();
47
53 void setup(IntVec<D> const & meshDimensions);
54
61 void setup(RDField<D>& rDField, RDFieldDft<D>& kDField);
62
69 void forwardTransform(RDField<D> & rField, RDFieldDft<D>& kField)
70 const;
71
78 void forwardTransformSafe(RDField<D> const & rField,
79 RDFieldDft<D>& kField) const;
80
87 void inverseTransform(RDFieldDft<D> & kField, RDField<D>& rField)
88 const;
89
96 void inverseTransformSafe(RDFieldDft<D> const & kField,
97 RDField<D>& rField) const;
98
102 const IntVec<D>& meshDimensions() const;
103
107 bool isSetup() const;
108
112 cufftHandle& fPlan();
113
117 cufftHandle& iPlan();
118
119 private:
120
121 // Vector containing number of grid points in each direction.
122 IntVec<D> meshDimensions_;
123
124 // Private r-space array for performing safe transforms.
125 mutable RDField<D> rFieldCopy_;
126
127 // Private k-space array for performing safe transforms.
128 mutable RDFieldDft<D> kFieldCopy_;
129
130 // Number of points in r-space grid
131 int rSize_;
132
133 // Number of points in k-space grid
134 int kSize_;
135
136 // Pointer to a plan for a forward transform.
137 cufftHandle fPlan_;
138
139 // Pointer to a plan for an inverse transform.
140 cufftHandle iPlan_;
141
142 // Have array dimension and plan been initialized?
143 bool isSetup_;
144
148 void makePlans(RDField<D>& rDField, RDFieldDft<D>& kDField);
149
150 };
151
152 // Declarations of explicit specializations
153
154 template <>
155 void FFT<1>::makePlans(RDField<1>& rField, RDFieldDft<1>& kField);
156
157 template <>
158 void FFT<2>::makePlans(RDField<2>& rField, RDFieldDft<2>& kField);
159
160 template <>
161 void FFT<3>::makePlans(RDField<3>& rField, RDFieldDft<3>& kField);
162
163 /*
164 * Return the dimensions of the grid for which this was allocated.
165 */
166 template <int D>
167 inline const IntVec<D>& FFT<D>::meshDimensions() const
168 { return meshDimensions_; }
169
170 template <int D>
171 inline bool FFT<D>::isSetup() const
172 { return isSetup_; }
173
174 template <int D>
175 inline cufftHandle& FFT<D>::fPlan()
176 { return fPlan_; }
177
178 template <int D>
179 inline cufftHandle& FFT<D>::iPlan()
180 { return iPlan_; }
181
182
183 #ifndef PSPG_FFT_TPP
184 // Suppress implicit instantiation
185 extern template class FFT<1>;
186 extern template class FFT<2>;
187 extern template class FFT<3>;
188 #endif
189
190static __global__
191void scaleRealData(cudaReal* data, cudaReal scale, int size) {
192 //write code that will scale
193 int nThreads = blockDim.x * gridDim.x;
194 int startId = blockIdx.x * blockDim.x + threadIdx.x;
195 for(int i = startId; i < size; i += nThreads ) {
196 data[i] *= scale;
197 }
198}
199
200}
201}
202
203//#include "FFT.tpp"
204#endif
An IntVec<D, T> is a D-component vector of elements of integer type T.
Definition: IntVec.h:27
Fourier transform wrapper for real data.
void forwardTransformSafe(RDField< D > const &rField, RDFieldDft< D > &kField) const
Compute forward Fourier transform without destroying input.
void forwardTransform(RDField< D > &rField, RDFieldDft< D > &kField) const
Compute forward (real-to-complex) discrete Fourier transform.
void setup(IntVec< D > const &meshDimensions)
Setup grid dimensions, plans and work space.
const IntVec< D > & meshDimensions() const
Return the dimensions of the grid for which this was allocated.
bool isSetup() const
Has this FFT object been setup?
cufftHandle & fPlan()
Get the plan for the forward DFT.
void inverseTransform(RDFieldDft< D > &kField, RDField< D > &rField) const
Compute inverse (complex-to-real) discrete Fourier transform.
cufftHandle & iPlan()
Get the plan for the inverse DFT.
virtual ~FFT()
Destructor.
FFT()
Default constructor.
void inverseTransformSafe(RDFieldDft< D > const &kField, RDField< D > &rField) const
Compute inverse (complex to real) DFT without destroying input.
Discrete Fourier Transform (DFT) of a real field on an FFT mesh.
Definition: RDFieldDft.h:35
Field of real single precision values on an FFT mesh on a device.
Definition: RDField.h:34
File containing preprocessor macros for error handling.
C++ namespace for polymer self-consistent field theory (PSCF).
Utility classes for scientific computation.
Definition: accumulators.mod:1