ls1-MarDyn
ls1-MarDyn molecular dynamics code
FFTW_API.h
1/*
2 * FFTW_API.h
3 *
4 * Created on: Sep 14, 2015
5 * Author: gallardjm
6 */
7
8#ifndef FFTW_API_H_
9#define FFTW_API_H_
10
11#include <stddef.h> //NULL
12#include <fftw3.h>
13#include "bhfmm/fft/FFTSettings_preprocessor.h"
14
25class FFTW_API {
26
27public:
28
35 FFTW_API(int nx, int ny) :
36 _nx(nx), _ny(ny) {
37#if defined(__SINGLE_PRECISION_FFT__)
38 _f_in = (fftwf_complex*)fftwf_malloc(sizeof(fftwf_complex)*nx*ny);
39 _f_out = (fftwf_complex*)fftwf_malloc(sizeof(fftwf_complex)*nx*ny);
40 _b_in = (fftwf_complex*)fftwf_malloc(sizeof(fftwf_complex)*nx*ny);
41 _b_out = (fftwf_complex*)fftwf_malloc(sizeof(fftwf_complex)*nx*ny);
42
43 _forward = fftwf_plan_dft_2d(nx, ny, _f_in, _f_out, FFTW_FORWARD, FFTW_ESTIMATE);
44 _backward = fftwf_plan_dft_2d(nx, ny, _b_in, _b_out, FFTW_BACKWARD, FFTW_ESTIMATE);
45#else
46 _f_in = (fftw_complex*) fftwf_malloc(sizeof(fftw_complex) * nx * ny);
47 _f_out = (fftw_complex*) fftwf_malloc(sizeof(fftw_complex) * nx * ny);
48 _b_in = (fftw_complex*) fftwf_malloc(sizeof(fftw_complex) * nx * ny);
49 _b_out = (fftw_complex*) fftwf_malloc(sizeof(fftw_complex) * nx * ny);
50
51 //if(ny>1) {
52 _forward = fftw_plan_dft_2d(nx, ny, _f_in, _f_out, FFTW_FORWARD,
53 FFTW_ESTIMATE);
54 _backward = fftw_plan_dft_2d(nx, ny, _b_in, _b_out, FFTW_BACKWARD,
55 FFTW_ESTIMATE);
56 //} else {
57 // _forward = fftw_plan_dft_1d(nx, _f_in, _f_out, FFTW_FORWARD, FFTW_ESTIMATE);
58 // _backward = fftw_plan_dft_1d(nx, _b_in, _b_out, FFTW_BACKWARD, FFTW_ESTIMATE);
59 //}
60#endif
61 }
62
65#if defined(__SINGLE_PRECISION_FFT__)
66 fftwf_destroy_plan(_forward);
67 fftwf_destroy_plan(_backward);
68 fftwf_free(_f_in);
69 fftwf_free(_f_out);
70 fftwf_free(_b_in);
71 fftwf_free(_b_out);
72#else
73 fftw_destroy_plan(_forward);
74 fftw_destroy_plan(_backward);
75 fftw_free(_f_in);
76 fftw_free(_f_out);
77 fftw_free(_b_in);
78 fftw_free(_b_out);
79
80 fftw_cleanup();
81#endif
82 }
83
84#if defined(__SINGLE_PRECISION_FFT__)
85 fftwf_complex* getIn_Forward() {return _f_in;}
86 fftwf_complex* getIn_Backward() {return _b_in;}
87
88 fftwf_complex* FFTAndGetOutput_Forward() {fftwf_execute(_forward); return _f_out;}
89 fftwf_complex* FFTAndGetOutput_Backward() {fftwf_execute(_backward); return _b_out;}
90#else
91 fftw_complex* getIn_Forward() {
92 return _f_in;
93 }
94 fftw_complex* getIn_Backward() {
95 return _b_in;
96 }
97
98 fftw_complex* FFTAndGetOutput_Forward() {
99 fftw_execute(_forward);
100 return _f_out;
101 }
102 fftw_complex* FFTAndGetOutput_Backward() {
103 fftw_execute(_backward);
104 return _b_out;
105 }
106#endif
107
108private:
109 int _nx;
110 int _ny;
111
112#if defined(__SINGLE_PRECISION_FFT__)
113 fftwf_plan _forward;
114 fftwf_complex* _f_in;
115 fftwf_complex* _f_out;
116 fftwf_plan _backward;
117 fftwf_complex* _b_in;
118 fftwf_complex* _b_out;
119#else
120 fftw_plan _forward;
121 fftw_complex* _f_in;
122 fftw_complex* _f_out;
123 fftw_plan _backward;
124 fftw_complex* _b_in;
125 fftw_complex* _b_out;
126#endif
127};
128
129#endif
Definition: FFTW_API.h:25
~FFTW_API()
destructor, clean plans and frees allocated memory
Definition: FFTW_API.h:64
FFTW_API(int nx, int ny)
Definition: FFTW_API.h:35