1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
/* -*- c++ -*- ----------------------------------------------------------
LAMMPS - Large-scale Atomic/Molecular Massively Parallel Simulator
http://lammps.sandia.gov, Sandia National Laboratories
Steve Plimpton, sjplimp@sandia.gov
Copyright (2003) Sandia Corporation. Under the terms of Contract
DE-AC04-94AL85000 with Sandia Corporation, the U.S. Government retains
certain rights in this software. This software is distributed under
the GNU General Public License.
See the README file in the top-level LAMMPS directory.
------------------------------------------------------------------------- */
#include <mpi.h>
// User-settable FFT precision
// FFT_PRECISION = 1 is single-precision complex (4-byte real, 4-byte imag)
// FFT_PRECISION = 2 is double-precision complex (8-byte real, 8-byte imag)
#ifdef FFT_SINGLE
#define FFT_PRECISION 1
typedef float FFT_SCALAR;
#else
#define FFT_PRECISION 2
typedef double FFT_SCALAR;
#endif
// if user set FFTW, it means FFTW3
#ifdef FFT_FFTW
#define FFT_FFTW3
#endif
// -------------------------------------------------------------------------
// Data types for single-precision complex
#if FFT_PRECISION == 1
#if defined(FFT_MKL)
#include "mkl_dfti.h"
typedef float _Complex FFT_DATA;
#define FFT_MKL_PREC DFTI_SINGLE
//#elif defined(FFT_FFTW2)
//#if defined(FFTW_SIZE)
//#include "sfftw.h"
//#else
//#include "fftw.h"
//#endif
//typedef FFTW_COMPLEX FFT_DATA;
#elif defined(FFT_FFTW3)
#include "fftw3.h"
typedef fftwf_complex FFT_DATA;
#define FFTW_API(function) fftwf_ ## function
#else
/* use a stripped down version of kiss fft as default fft */
#ifndef FFT_KISS
#define FFT_KISS
#endif
#define kiss_fft_scalar float
typedef struct {
kiss_fft_scalar re;
kiss_fft_scalar im;
} FFT_DATA;
struct kiss_fft_state;
typedef struct kiss_fft_state* kiss_fft_cfg;
#endif
// -------------------------------------------------------------------------
// Data types for double-precision complex
#elif FFT_PRECISION == 2
#if defined(FFT_MKL)
#include "mkl_dfti.h"
typedef double _Complex FFT_DATA;
#define FFT_MKL_PREC DFTI_DOUBLE
//#elif defined(FFT_FFTW2)
//#if defined(FFTW_SIZE)
//#include "dfftw.h"
//#else
//#include "fftw.h"
//#endif
//typedef FFTW_COMPLEX FFT_DATA;
#elif defined(FFT_FFTW3)
#include "fftw3.h"
typedef fftw_complex FFT_DATA;
#define FFTW_API(function) fftw_ ## function
#else
/* use a stripped down version of kiss fft as default fft */
#ifndef FFT_KISS
#define FFT_KISS
#endif
#define kiss_fft_scalar double
typedef struct {
kiss_fft_scalar re;
kiss_fft_scalar im;
} FFT_DATA;
struct kiss_fft_state;
typedef struct kiss_fft_state* kiss_fft_cfg;
#endif
#else
#error "FFT_PRECISION needs to be either 1 (=single) or 2 (=double)"
#endif
// -------------------------------------------------------------------------
// details of how to do a 3d FFT
struct fft_plan_3d {
struct remap_plan_3d *pre_plan; // remap from input -> 1st FFTs
struct remap_plan_3d *mid1_plan; // remap from 1st -> 2nd FFTs
struct remap_plan_3d *mid2_plan; // remap from 2nd -> 3rd FFTs
struct remap_plan_3d *post_plan; // remap from 3rd FFTs -> output
FFT_DATA *copy; // memory for remap results (if needed)
FFT_DATA *scratch; // scratch space for remaps
int total1,total2,total3; // # of 1st,2nd,3rd FFTs (times length)
int length1,length2,length3; // length of 1st,2nd,3rd FFTs
int pre_target; // where to put remap results
int mid1_target,mid2_target;
int scaled; // whether to scale FFT results
int normnum; // # of values to rescale
double norm; // normalization factor for rescaling
// system specific 1d FFT info
#if defined(FFT_MKL)
DFTI_DESCRIPTOR *handle_fast;
DFTI_DESCRIPTOR *handle_mid;
DFTI_DESCRIPTOR *handle_slow;
//#elif defined(FFT_FFTW2)
// fftw_plan plan_fast_forward;
// fftw_plan plan_fast_backward;
// fftw_plan plan_mid_forward;
// fftw_plan plan_mid_backward;
//fftw_plan plan_slow_forward;
//fftw_plan plan_slow_backward;
#elif defined(FFT_FFTW3)
FFTW_API(plan) plan_fast_forward;
FFTW_API(plan) plan_fast_backward;
FFTW_API(plan) plan_mid_forward;
FFTW_API(plan) plan_mid_backward;
FFTW_API(plan) plan_slow_forward;
FFTW_API(plan) plan_slow_backward;
#elif defined(FFT_KISS)
kiss_fft_cfg cfg_fast_forward;
kiss_fft_cfg cfg_fast_backward;
kiss_fft_cfg cfg_mid_forward;
kiss_fft_cfg cfg_mid_backward;
kiss_fft_cfg cfg_slow_forward;
kiss_fft_cfg cfg_slow_backward;
#endif
};
// function prototypes
extern "C" {
void fft_3d(FFT_DATA *, FFT_DATA *, int, struct fft_plan_3d *);
struct fft_plan_3d *fft_3d_create_plan(MPI_Comm, int, int, int,
int, int, int, int, int,
int, int, int, int, int, int, int,
int, int, int *, int);
void fft_3d_destroy_plan(struct fft_plan_3d *);
void factor(int, int *, int *);
void bifactor(int, int *, int *);
void fft_1d_only(FFT_DATA *, int, int, struct fft_plan_3d *);
}
/* ERROR/WARNING messages:
*/