PSCF v1.1
FFTBatched.h
1#ifndef PSPG_FFT_BATCHED_H
2#define PSPG_FFT_BATCHED_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 <pscf/math/IntVec.h>
14#include <util/global.h>
15
16#include <cufft.h>
17#include <cuda.h>
18#include <cuda_runtime.h>
19#include <pspg/math/GpuResources.h>
20
21//temporary for debugging
22#include <iostream>
23
24namespace Pscf {
25namespace Pspg {
26
27 using namespace Util;
28 using namespace Pscf;
29
35 template <int D>
37 {
38
39 public:
40
44 FFTBatched();
45
49 virtual ~FFTBatched();
50
57 void setup(RDField<D>& rDField, RDFieldDft<D>& kDField);
58
59 void setup(const IntVec<D>& rDim,const IntVec<D>& kDim, int batchSize);
67
68 void forwardTransform(cudaReal* in, cudaComplex* out, int batchSize);
76
77 void inverseTransform(cudaComplex* in, cudaReal* out, int batchSize);
78
82 const IntVec<D>& meshDimensions() const;
83
84 bool isSetup() const;
85
86 cufftHandle& fPlan();
87
88 cufftHandle& iPlan();
89
90
91 private:
92
93 // Vector containing number of grid points in each direction.
94 IntVec<D> meshDimensions_;
95
96 // Number of points in r-space grid
97 int rSize_;
98
99 // Number of points in k-space grid
100 int kSize_;
101
102 // Pointer to a plan for a forward transform.
103 cufftHandle fPlan_;
104
105 // Pointer to a plan for an inverse transform.
106 cufftHandle iPlan_;
107
108 // Have array dimension and plan been initialized?
109 bool isSetup_;
110
114 void makePlans(RDField<D>& rDField, RDFieldDft<D>& kDField);
115 void makePlans(const IntVec<D>& rDim, const IntVec<D>& kDField, int batchSize);
116
117 };
118
119 /*
120 * Return the dimensions of the grid for which this was allocated.
121 */
122 template <int D>
124 { return meshDimensions_; }
125
126 template <int D>
127 inline bool FFTBatched<D>::isSetup() const
128 { return isSetup_; }
129
130 template <int D>
131 inline cufftHandle& FFTBatched<D>::fPlan()
132 { return fPlan_; }
133
134 template <int D>
135 inline cufftHandle& FFTBatched<D>::iPlan()
136 { return iPlan_; }
137
138 #ifndef PSPG_FFT_BATCHED_TPP
139 // Suppress implicit instantiation
140 extern template class FFTBatched<1>;
141 extern template class FFTBatched<2>;
142 extern template class FFTBatched<3>;
143 #endif
144}
145}
146//#include "FFTBatched.tpp"
147#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.
Definition: FFTBatched.h:37
FFTBatched()
Default constructor.
Definition: FFTBatched.tpp:46
void setup(RDField< D > &rDField, RDFieldDft< D > &kDField)
Check and setup grid dimensions if necessary.
Definition: FFTBatched.tpp:73
void inverseTransform(RDFieldDft< D > &in, RDField< D > &out)
Compute inverse (complex-to-real) Fourier transform.
Definition: FFTBatched.tpp:326
virtual ~FFTBatched()
Destructor.
Definition: FFTBatched.tpp:59
void forwardTransform(RDField< D > &in, RDFieldDft< D > &out)
Compute forward (real-to-complex) Fourier transform.
Definition: FFTBatched.tpp:248
const IntVec< D > & meshDimensions() const
Return the dimensions of the grid for which this was allocated.
Definition: FFTBatched.h:123
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