Bridge++  Ver. 2.0.2
fft_xyz_3dim.cpp
Go to the documentation of this file.
1 
14 #ifdef USE_FFTWLIB
15 
16 #include "fft_xyz_3dim.h"
17 
18 #ifdef USE_FACTORY_AUTOREGISTER
19 namespace {
20  bool init = FFT_xyz_3dim::register_factory();
21 }
22 #endif
23 
24 const std::string FFT_xyz_3dim::class_name = "FFT_xyz_3dim";
25 
26 //====================================================================
27 void FFT_xyz_3dim::set_parameters(const Parameters& params)
28 {
29  std::string vlevel;
30  if (!params.fetch_string("verbose_level", vlevel)) {
31  m_vl = vout.set_verbose_level(vlevel);
32  }
33 
34  //- fetch and check input parameters
35  string str_fft_direction;
36 
37  int err = 0;
38  err += params.fetch_string("FFT_direction", str_fft_direction);
39 
40  if (err) {
41  vout.crucial(m_vl, "Error at %s: input parameter not found.\n", class_name.c_str());
42  exit(EXIT_FAILURE);
43  }
44 
45  set_parameters(str_fft_direction);
46 }
47 
48 
49 //====================================================================
50 void FFT_xyz_3dim::get_parameters(Parameters& params) const
51 {
52  params.set_string("FFT_direction", m_is_forward ? "Forward" : "Backward");
53 
54  params.set_string("verbose_level", vout.get_verbose_level(m_vl));
55 }
56 
57 
58 //====================================================================
59 void FFT_xyz_3dim::set_parameters(const string& str_fft_direction)
60 {
61  //- print input parameters
62  vout.general(m_vl, "%s:\n", class_name.c_str());
63  vout.general(m_vl, " FFT_direction = %s\n", str_fft_direction.c_str());
64 
65  //- range check
66 
67  //- store values
68  if (str_fft_direction == "Forward") {
69  m_is_forward = true;
70  } else if (str_fft_direction == "Backward") {
71  m_is_forward = false;
72  } else {
73  vout.crucial(m_vl, "Error at %s: unsupported FFT direction \"%s\"\n", class_name.c_str(), str_fft_direction.c_str());
74  exit(EXIT_FAILURE);
75  }
76 }
77 
78 
79 //====================================================================
80 void FFT_xyz_3dim::init()
81 {
82  //- global lattice size
83  const int Lx = CommonParameters::Lx();
84  const int Ly = CommonParameters::Ly();
85  const int Lz = CommonParameters::Lz();
86 
87 #ifdef USE_OPENMP
88  int threads_ok = fftw_init_threads();
89 #endif
90 
91 #ifdef USE_MPI
92  const int NPE_x = CommonParameters::NPEx();
93  const int NPE_y = CommonParameters::NPEy();
94  // const int NPE_z = CommonParameters::NPEz();
95  const int NPE_t = CommonParameters::NPEt();
96 
97  if ((NPE_x * NPE_y * NPE_t) != 1) {
98  vout.crucial(m_vl, "Error at %s: FFTW supports parallelization only in z-direction.\n",
99  class_name.c_str());
100  exit(EXIT_FAILURE);
101  }
102 
103 
104  fftw_mpi_init();
105 
106 
107  //- allocate m_in,out = m_in,out[Nz][Ly][Lx]
108  const ptrdiff_t Lx_p = CommonParameters::Lx();
109  const ptrdiff_t Ly_p = CommonParameters::Ly();
110  const ptrdiff_t Lz_p = CommonParameters::Lz();
111 
112  ptrdiff_t fftw_size_p = fftw_mpi_local_size_3d(Lz_p, Ly_p, Lx_p,
113  Communicator_impl::world(),
114  &m_Nz_p, &m_z_start_p);
115 
116  m_in = fftw_alloc_complex(fftw_size_p);
117  m_out = fftw_alloc_complex(fftw_size_p);
118 
119  if (!m_in || !m_out) {
120  vout.crucial(m_vl, "Error at %s: failed to allocate memory %d [Byte].\n",
121  class_name.c_str(), (int)fftw_size_p);
122  exit(EXIT_FAILURE);
123  }
124 #else
125  //- allocate m_in,out = m_in,out[Nz][Ly][Lx]
126  const size_t fftw_size = sizeof(fftw_complex) * Lx * Ly * Lz;
127  m_in = (fftw_complex *)fftw_malloc(fftw_size);
128  m_out = (fftw_complex *)fftw_malloc(fftw_size);
129 
130  if (!m_in || !m_out) {
131  vout.crucial(m_vl, "Error at %s: failed to allocate memory %d [Byte].\n",
132  class_name.c_str(), (int)fftw_size);
133  exit(EXIT_FAILURE);
134  }
135 #endif
136 }
137 
138 
139 //====================================================================
140 void FFT_xyz_3dim::tidy_up()
141 {
142  if (m_in) fftw_free(m_in);
143  if (m_out) fftw_free(m_out);
144  if (m_plan) fftw_destroy_plan(m_plan);
145 }
146 
147 
148 //====================================================================
149 void FFT_xyz_3dim::fft(Field& field)
150 {
151  //- global lattice size
152  const int Lx = CommonParameters::Lx();
153  const int Ly = CommonParameters::Ly();
154  const int Lz = CommonParameters::Lz();
155  const int Lt = CommonParameters::Lt();
156  const int Lxyz = Lx * Ly * Lz;
157 
158  //- local size
159  const int Nz = CommonParameters::Nz();
160 
161  const int Nin = field.nin();
162  const int Nex = field.nex();
163 
164 
165  //- setup FFTW plan
166 #ifdef USE_OPENMP
167  const int Nthread = ThreadManager::get_num_threads();
168  fftw_plan_with_nthreads(Nthread);
169 #endif
170 #ifdef USE_MPI
171  const ptrdiff_t Lx_p = CommonParameters::Lx();
172  const ptrdiff_t Ly_p = CommonParameters::Ly();
173  const ptrdiff_t Lz_p = CommonParameters::Lz();
174 
175  if (m_is_forward) {
176  m_plan = fftw_mpi_plan_dft_3d(Lz_p, Ly_p, Lx_p, m_in, m_out,
177  Communicator_impl::world(),
178  FFTW_FORWARD, FFTW_ESTIMATE);
179  } else {
180  m_plan = fftw_mpi_plan_dft_3d(Lz_p, Ly_p, Lx_p, m_in, m_out,
181  Communicator_impl::world(),
182  FFTW_BACKWARD, FFTW_ESTIMATE);
183  }
184 #else
185  if (m_is_forward) {
186  m_plan = fftw_plan_dft_3d(Lz, Ly, Lx, m_in, m_out,
187  FFTW_FORWARD, FFTW_ESTIMATE);
188  } else {
189  m_plan = fftw_plan_dft_3d(Lz, Ly, Lx, m_in, m_out,
190  FFTW_BACKWARD, FFTW_ESTIMATE);
191  }
192 #endif
193 
194 
195  // #### Execution main part ####
196  //- Nin is devided by 2, because of complex(i.e. real and imag)
197  for (int in2 = 0; in2 < Nin / 2; ++in2) {
198  for (int t_global = 0; t_global < Lt; t_global++) {
199  for (int ex = 0; ex < Nex; ++ex) {
200  //- input data
201  for (int z = 0; z < Nz; z++) {
202  for (int y_global = 0; y_global < Ly; y_global++) {
203  for (int x_global = 0; x_global < Lx; x_global++) {
204  int isite_xyz_local = x_global + Lx * (y_global + Ly * z);
205 
206  int isite = m_index.site(x_global, y_global, z, t_global);
207  int i_real = 2 * in2;
208  int i_imag = 2 * in2 + 1;
209 
210  m_in[isite_xyz_local][0] = field.cmp(i_real, isite, ex);
211  m_in[isite_xyz_local][1] = field.cmp(i_imag, isite, ex);
212  }
213  }
214  }
215 
216 
217  fftw_execute(m_plan);
218 
219 
220  //- output data
221  for (int z = 0; z < Nz; z++) {
222  for (int y_global = 0; y_global < Ly; y_global++) {
223  for (int x_global = 0; x_global < Lx; x_global++) {
224  int isite_xyz_local = x_global + Lx * (y_global + Ly * z);
225 
226  int isite = m_index.site(x_global, y_global, z, t_global);
227  int i_real = 2 * in2;
228  int i_imag = 2 * in2 + 1;
229 
230  field.set(i_real, isite, ex, m_out[isite_xyz_local][0]);
231  field.set(i_imag, isite, ex, m_out[isite_xyz_local][1]);
232  }
233  }
234  }
235  }
236  }
237  }
238  //- end of global loops
239 
240  //- normailzation for FFTW_BACKWARD
241  if (!m_is_forward) {
242  scal(field, 1.0 / Lxyz);
243  }
244 }
245 
246 
247 //====================================================================
248 void FFT_xyz_3dim::fft(Field& field_out, const Field& field_in)
249 {
250  //- global lattice size
251  const int Lx = CommonParameters::Lx();
252  const int Ly = CommonParameters::Ly();
253  const int Lz = CommonParameters::Lz();
254  const int Lt = CommonParameters::Lt();
255  const int Lxyz = Lx * Ly * Lz;
256 
257  //- local size
258  const int Nz = CommonParameters::Nz();
259 
260  const int Nin = field_in.nin();
261  const int Nex = field_in.nex();
262 
263 
264  //- setup FFTW plan
265 #ifdef USE_OPENMP
266  const int Nthread = ThreadManager::get_num_threads();
267  fftw_plan_with_nthreads(Nthread);
268 #endif
269 #ifdef USE_MPI
270  const ptrdiff_t Lx_p = CommonParameters::Lx();
271  const ptrdiff_t Ly_p = CommonParameters::Ly();
272  const ptrdiff_t Lz_p = CommonParameters::Lz();
273 
274  if (m_is_forward) {
275  m_plan = fftw_mpi_plan_dft_3d(Lz_p, Ly_p, Lx_p, m_in, m_out,
276  Communicator_impl::world(),
277  FFTW_FORWARD, FFTW_ESTIMATE);
278  } else {
279  m_plan = fftw_mpi_plan_dft_3d(Lz_p, Ly_p, Lx_p, m_in, m_out,
280  Communicator_impl::world(),
281  FFTW_BACKWARD, FFTW_ESTIMATE);
282  }
283 #else
284  if (m_is_forward) {
285  m_plan = fftw_plan_dft_3d(Lz, Ly, Lx, m_in, m_out,
286  FFTW_FORWARD, FFTW_ESTIMATE);
287  } else {
288  m_plan = fftw_plan_dft_3d(Lz, Ly, Lx, m_in, m_out,
289  FFTW_BACKWARD, FFTW_ESTIMATE);
290  }
291 #endif
292 
293 
294  // #### Execution main part ####
295  //- Nin is devided by 2, because of complex(i.e. real and imag)
296  for (int in2 = 0; in2 < Nin / 2; ++in2) {
297  for (int t_global = 0; t_global < Lt; t_global++) {
298  for (int ex = 0; ex < Nex; ++ex) {
299  //- input data
300  for (int z = 0; z < Nz; z++) {
301  for (int y_global = 0; y_global < Ly; y_global++) {
302  for (int x_global = 0; x_global < Lx; x_global++) {
303  int isite_xyz_local = x_global + Lx * (y_global + Ly * z);
304 
305  int isite = m_index.site(x_global, y_global, z, t_global);
306  int i_real = 2 * in2;
307  int i_imag = 2 * in2 + 1;
308 
309  m_in[isite_xyz_local][0] = field_in.cmp(i_real, isite, ex);
310  m_in[isite_xyz_local][1] = field_in.cmp(i_imag, isite, ex);
311  }
312  }
313  }
314 
315 
316  fftw_execute(m_plan);
317 
318 
319  //- output data
320  for (int z = 0; z < Nz; z++) {
321  for (int y_global = 0; y_global < Ly; y_global++) {
322  for (int x_global = 0; x_global < Lx; x_global++) {
323  int isite_xyz_local = x_global + Lx * (y_global + Ly * z);
324 
325  int isite = m_index.site(x_global, y_global, z, t_global);
326  int i_real = 2 * in2;
327  int i_imag = 2 * in2 + 1;
328 
329  field_out.set(i_real, isite, ex, m_out[isite_xyz_local][0]);
330  field_out.set(i_imag, isite, ex, m_out[isite_xyz_local][1]);
331  }
332  }
333  }
334  }
335  }
336  }
337  //- end of global loops
338 
339  //- normailzation for FFTW_BACKWARD
340  if (!m_is_forward) {
341  scal(field_out, 1.0 / Lxyz);
342  }
343 }
344 
345 
346 //====================================================================
347 void FFT_xyz_3dim::fft(Field& field_out, const Field& field_in, const Direction dir)
348 {
349  // save state
350  bool backup_fwbw = m_is_forward;
351 
352  // find direction and set
353  if (dir == FORWARD) {
354  m_is_forward = true;
355  } else if (dir == BACKWARD) {
356  m_is_forward = false;
357  } else {
358  vout.crucial(m_vl, "%s: unknown direction %d. failed.\n", class_name.c_str(), dir);
359  exit(EXIT_FAILURE);
360  }
361 
362  // delegate to another method
363  fft(field_out, field_in);
364 
365  // restore state
366  m_is_forward = backup_fwbw;
367 }
368 
369 
370 //==========================================================
371 //==================================================END=====
372 #endif
CommonParameters::Nz
static int Nz()
Definition: commonParameters.h:107
Parameters::set_string
void set_string(const string &key, const string &value)
Definition: parameters.cpp:39
ThreadManager::get_num_threads
static int get_num_threads()
returns available number of threads.
Definition: threadManager.cpp:246
Field::set
void set(const int jin, const int site, const int jex, double v)
Definition: field.h:175
Parameters
Class for parameters.
Definition: parameters.h:46
Field::nex
int nex() const
Definition: field.h:128
CommonParameters::Ly
static int Ly()
Definition: commonParameters.h:92
Direction
Direction
Definition: bridge_defs.h:24
Field::nin
int nin() const
Definition: field.h:126
CommonParameters::Lt
static int Lt()
Definition: commonParameters.h:94
CommonParameters::Lx
static int Lx()
Definition: commonParameters.h:91
CommonParameters::Lz
static int Lz()
Definition: commonParameters.h:93
CommonParameters::NPEy
static int NPEy()
Definition: commonParameters.h:98
Field::cmp
double cmp(const int jin, const int site, const int jex) const
Definition: field.h:143
Bridge::BridgeIO::set_verbose_level
static VerboseLevel set_verbose_level(const std::string &str)
Definition: bridgeIO.cpp:133
CommonParameters::NPEx
static int NPEx()
Definition: commonParameters.h:97
scal
void scal(Field &x, const double a)
scal(x, a): x = a * x
Definition: field.cpp:261
Parameters::fetch_string
int fetch_string(const string &key, string &value) const
Definition: parameters.cpp:378
CommonParameters::NPEt
static int NPEt()
Definition: commonParameters.h:100
Bridge::BridgeIO::crucial
void crucial(const char *format,...)
Definition: bridgeIO.cpp:180
Field
Container of Field-type object.
Definition: field.h:46
fft_xyz_3dim.h
Bridge::BridgeIO::general
void general(const char *format,...)
Definition: bridgeIO.cpp:200
Bridge::vout
BridgeIO vout
Definition: bridgeIO.cpp:512
Bridge::BridgeIO::get_verbose_level
static std::string get_verbose_level(const VerboseLevel vl)
Definition: bridgeIO.cpp:154