18 #ifdef USE_FACTORY_AUTOREGISTER
20 bool init = FFT_xyz_3dim::register_factory();
24 const std::string FFT_xyz_3dim::class_name =
"FFT_xyz_3dim";
27 void FFT_xyz_3dim::set_parameters(
const Parameters& params)
35 string str_fft_direction;
38 err += params.
fetch_string(
"FFT_direction", str_fft_direction);
41 vout.
crucial(m_vl,
"Error at %s: input parameter not found.\n", class_name.c_str());
45 set_parameters(str_fft_direction);
50 void FFT_xyz_3dim::get_parameters(
Parameters& params)
const
52 params.
set_string(
"FFT_direction", m_is_forward ?
"Forward" :
"Backward");
59 void FFT_xyz_3dim::set_parameters(
const string& str_fft_direction)
63 vout.
general(m_vl,
" FFT_direction = %s\n", str_fft_direction.c_str());
68 if (str_fft_direction ==
"Forward") {
70 }
else if (str_fft_direction ==
"Backward") {
73 vout.
crucial(m_vl,
"Error at %s: unsupported FFT direction \"%s\"\n", class_name.c_str(), str_fft_direction.c_str());
80 void FFT_xyz_3dim::init()
88 int threads_ok = fftw_init_threads();
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",
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);
116 m_in = fftw_alloc_complex(fftw_size_p);
117 m_out = fftw_alloc_complex(fftw_size_p);
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);
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);
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);
140 void FFT_xyz_3dim::tidy_up()
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);
149 void FFT_xyz_3dim::fft(
Field& field)
156 const int Lxyz = Lx * Ly * Lz;
161 const int Nin = field.
nin();
162 const int Nex = field.
nex();
168 fftw_plan_with_nthreads(Nthread);
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);
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);
186 m_plan = fftw_plan_dft_3d(Lz, Ly, Lx, m_in, m_out,
187 FFTW_FORWARD, FFTW_ESTIMATE);
189 m_plan = fftw_plan_dft_3d(Lz, Ly, Lx, m_in, m_out,
190 FFTW_BACKWARD, FFTW_ESTIMATE);
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) {
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);
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;
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);
217 fftw_execute(m_plan);
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);
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;
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]);
242 scal(field, 1.0 / Lxyz);
248 void FFT_xyz_3dim::fft(
Field& field_out,
const Field& field_in)
255 const int Lxyz = Lx * Ly * Lz;
260 const int Nin = field_in.
nin();
261 const int Nex = field_in.
nex();
267 fftw_plan_with_nthreads(Nthread);
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);
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);
285 m_plan = fftw_plan_dft_3d(Lz, Ly, Lx, m_in, m_out,
286 FFTW_FORWARD, FFTW_ESTIMATE);
288 m_plan = fftw_plan_dft_3d(Lz, Ly, Lx, m_in, m_out,
289 FFTW_BACKWARD, FFTW_ESTIMATE);
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) {
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);
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;
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);
316 fftw_execute(m_plan);
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);
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;
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]);
341 scal(field_out, 1.0 / Lxyz);
350 bool backup_fwbw = m_is_forward;
353 if (dir == FORWARD) {
355 }
else if (dir == BACKWARD) {
356 m_is_forward =
false;
358 vout.
crucial(m_vl,
"%s: unknown direction %d. failed.\n", class_name.c_str(), dir);
363 fft(field_out, field_in);
366 m_is_forward = backup_fwbw;