OPAL (Object Oriented Parallel Accelerator Library)  2.2.0
OPAL
FFTBase.h
Go to the documentation of this file.
1 //
2 // IPPL FFT
3 //
4 // Copyright (c) 2008-2018
5 // Paul Scherrer Institut, Villigen PSI, Switzerland
6 // All rights reserved.
7 //
8 // OPAL is licensed under GNU GPL version 3.
9 //
10 
11 //--------------------------------------------------------------------------
12 // Class FFTBase
13 //--------------------------------------------------------------------------
14 
15 #ifndef IPPL_FFT_FFTBASE_H
16 #define IPPL_FFT_FFTBASE_H
17 
18 // include files
19 #include "Utility/PAssert.h"
20 #include "Index/NDIndex.h"
21 #include "Field/GuardCellSizes.h"
22 
23 #include "FFT/fftpack_FFT.h"
24 
25 #include <map>
26 #include <iostream>
27 
28 // forward declarations
29 template <unsigned Dim, class T> class FFTBase;
30 template <unsigned Dim, class T>
31 std::ostream& operator<<(std::ostream&, const FFTBase<Dim,T>&);
32 
34 inline
35 std::string getTransformType(unsigned int i)
36 {
37  static const char* transformTypeString_g[4] = { "complex-to-complex FFT",
38  "real-to-complex FFT",
39  "sine transform",
40  "cosine transform" };
41 
42  return std::string(transformTypeString_g[i % 4]);
43 }
44 
53 template <unsigned Dim, class T>
54 class FFTBase {
55 
56 public:
57  // Some externally visible typedefs and enums.
58  enum { dimensions = Dim }; // dimension
59  typedef T Precision_t; // precision
60  typedef NDIndex<Dim> Domain_t; // domain type
61 
62  // Enumeration of transform types, used by derived FFT classes
64 
65  // Type used for performing 1D FFTs
67 
68  FFTBase() {}
69 
81  FFTBase(FFT_e transform, const Domain_t& domain,
82  const bool transformTheseDims[Dim], bool compressTemps);
83 
92  FFTBase(FFT_e transform, const Domain_t& domain, bool compressTemps);
93 
94  // destructor
95  virtual ~FFTBase(void) { delete [] activeDims_m; }
96 
102  void write(std::ostream& out) const;
103 
110  void setDirectionName(int direction, const char* directionName);
111 
117  void setNormFact(Precision_t nf) { normFact_m = nf; }
118 
125  int transVnodes() const {
127  return Ippl::maxFFTNodes();
128  else
129  return (-1);
130  }
131 
132 protected:
133 
140 
142  int getDirection(const char* directionName) const;
143 
145  bool transformDim(unsigned d) const;
146 
148  unsigned numTransformDims(void) const { return nTransformDims_m; }
149 
151  unsigned activeDimension(unsigned d) const;
152 
155 
157  Precision_t& getNormFact(void) { return normFact_m; }
158 
160  const Domain_t& getDomain(void) const { return Domain_m; }
161 
163  bool checkDomain(const Domain_t& dom1, const Domain_t& dom2) const;
164 
166  bool compressTemps(void) const { return compressTempFields_m; }
167 
168 private:
169 
171  std::map<const char*,int> directions_m;
172 
175  unsigned nTransformDims_m;
176  unsigned* activeDims_m;
177 
180 
183 
186 
189 };
190 
191 
192 // Inline function definitions
193 
195 template <unsigned Dim, class T>
196 inline std::ostream&
197 operator<<(std::ostream& out, const FFTBase<Dim,T>& fft)
198 {
199  fft.write(out);
200  return out;
201 }
202 
207 template <unsigned Dim, class T>
208 inline void
210  const char* directionName) {
211  PAssert_EQ(std::abs(direction), 1);
212  directions_m[directionName] = direction;
213  return;
214 }
215 
223 template <unsigned Dim, class T>
224 inline int
225 FFTBase<Dim,T>::getDirection(const char* directionName) const {
226  return (*(directions_m.find(directionName))).second;
227 }
228 
236 template <unsigned Dim, class T>
237 inline bool
238 FFTBase<Dim,T>::transformDim(unsigned d) const {
239  PAssert_LT(d, Dim);
240  return transformDims_m[d];
241 }
242 
250 template <unsigned Dim, class T>
251 inline unsigned
253  PAssert_LT(d, nTransformDims_m);
254  return activeDims_m[d];
255 }
256 
267 template <unsigned Dim, class T>
268 inline bool
270  const FFTBase<Dim,T>::Domain_t& dom2) const {
271  // check whether domains are equivalent
272  // we require that some permutation of the axes gives a matching domain.
273  static bool matched[Dim];
274  bool found;
275  unsigned d, d1;
276  // initialize matched array to false
277  for (d=0; d<Dim; ++d) matched[d] = false;
278  d=0;
279  while (d<Dim) {
280  d1=0;
281  found = false;
282  while (!found && d1<Dim) {
283  // if we have not yet found a match for this dimension,
284  // compare length and base of Index objects
285  if (!matched[d1]) {
286  found = ( dom1[d].length()==dom2[d1].length() &&
287  dom1[d].sameBase(dom2[d1]) );
288  // if equivalent, mark this dimension as matched
289  if (found) matched[d1] = true;
290  }
291  ++d1;
292  }
293  if (!found) return false;
294  ++d;
295  }
296  return true;
297 }
298 
299 #include "FFT/FFTBase.hpp"
300 
301 #endif // IPPL_FFT_FFTBASE_H
302 
303 // vi: set et ts=4 sw=4 sts=4:
304 // Local Variables:
305 // mode:c
306 // c-basic-offset: 4
307 // indent-tabs-mode:nil
308 // End:
static int getNodes()
Definition: IpplInfo.cpp:773
FFTPACK< T > InternalFFT_t
Definition: FFTBase.h:66
Precision_t normFact_m
Normalization factor:
Definition: FFTBase.h:182
PETE_TUTree< FnAbs, typename T::PETE_Expr_t > abs(const PETE_Expr< T > &l)
std::map< const char *, int > directions_m
Stores user-defined names for FFT directions:
Definition: FFTBase.h:171
Definition: rbendmap.h:8
unsigned activeDimension(unsigned d) const
get dimension number from list of transformed dimensions
Definition: FFTBase.h:252
Domain_t Domain_m
Domain of the input field, mainly used to check axis sizes and ordering, former const Domain_t&amp; Domai...
Definition: FFTBase.h:185
const Domain_t & getDomain(void) const
get our domain
Definition: FFTBase.h:160
bool transformDim(unsigned d) const
query whether this dimension is to be transformed
Definition: FFTBase.h:238
bool compressTempFields_m
Switch to turn on/off compression of intermediate Fields (tempFields) as algorithm is finished with t...
Definition: FFTBase.h:188
NDIndex< Dim > Domain_t
Definition: FFTBase.h:60
FFTBase()
Definition: FFTBase.h:68
int transVnodes() const
Definition: FFTBase.h:125
InternalFFT_t FFTEngine_m
Internal FFT object for performing serial FFTs.
Definition: FFTBase.h:179
T Precision_t
Definition: FFTBase.h:59
#define PAssert_LT(a, b)
Definition: PAssert.h:121
unsigned numTransformDims(void) const
query number of transform dimensions
Definition: FFTBase.h:148
#define PAssert_EQ(a, b)
Definition: PAssert.h:119
FFT_e transformType_m
Indicates which type of transform we do.
Definition: FFTBase.h:173
bool transformDims_m[Dim]
Indicates which dimensions are transformed.
Definition: FFTBase.h:174
std::string getTransformType(unsigned int i)
character strings for transform types
Definition: FFTBase.h:35
InternalFFT_t & getEngine(void)
access the internal FFT Engine
Definition: FFTBase.h:154
virtual ~FFTBase(void)
Definition: FFTBase.h:95
static int maxFFTNodes()
Definition: IpplInfo.h:232
Precision_t & getNormFact(void)
get the FFT normalization factor
Definition: FFTBase.h:157
static GuardCellSizes< Dim > nullGC
null GuardCellSizes object for checking BareField arguments to transform
Definition: FFTBase.h:139
void setDirectionName(int direction, const char *directionName)
Definition: FFTBase.h:209
unsigned nTransformDims_m
Stores the number of dims to be transformed.
Definition: FFTBase.h:175
bool compressTemps(void) const
do we compress temps?
Definition: FFTBase.h:166
const unsigned Dim
bool checkDomain(const Domain_t &dom1, const Domain_t &dom2) const
compare indexes of two domains
Definition: FFTBase.h:269
unsigned * activeDims_m
Stores the numbers of these dims (0,1,2).
Definition: FFTBase.h:176
int getDirection(const char *directionName) const
translate direction name string into dimension number
Definition: FFTBase.h:225
void write(std::ostream &out) const
Definition: FFTBase.hpp:82
void setNormFact(Precision_t nf)
Definition: FFTBase.h:117