Fixed include statements such that double precision version of genborn.c
[gromacs/rigid-bodies.git] / src / mdlib / fft5d.h
blobd4444f268bee849c4d3c64449f2da1b9487066ce
1 #ifndef FFT5D_H_
2 #define FFT5D_H_
4 #ifdef HAVE_CONFIG_H
5 #include <config.h>
6 #endif
8 #ifdef NOGMX
9 /*#define GMX_MPI*/
10 /*#define GMX_FFT_FFTW3*/
11 FILE* debug;
12 #endif
14 #include <types/commrec.h>
15 #include "gmxcomplex.h"
16 #include "gmx_fft.h"
18 #ifndef GMX_LIB_MPI
19 double MPI_Wtime();
20 #endif
22 /*currently only special optimization for FFTE*/
23 #ifdef GMX_FFT_FFTW3
24 #include <fftw3.h>
25 #endif
27 #ifndef GMX_DOUBLE
28 #define FFTW(x) fftwf_##x
29 #else
30 #define FFTW(x) fftw_##x
31 #endif
33 struct fft5d_time_t {
34 double fft,local,mpi1,mpi2;
36 typedef struct fft5d_time_t *fft5d_time;
38 typedef enum fft5d_flags_t {
39 FFT5D_ORDER_YZ=1,
40 FFT5D_BACKWARD=2,
41 FFT5D_REALCOMPLEX=4,
42 FFT5D_DEBUG=8,
43 FFT5D_NOMEASURE=16,
44 FFT5D_INPLACE=32,
45 FFT5D_NOMALLOC=64
46 } fft5d_flags;
48 struct fft5d_plan_t {
49 t_complex *lin;
50 t_complex *lout;
51 gmx_fft_t p1d[3]; /*1D plans*/
52 #ifdef GMX_FFT_FFTW3
53 FFTW(plan) p2d; /*2D plan: used for 1D decomposition if FFT supports transposed output*/
54 FFTW(plan) p3d; /*3D plan: used for 0D decomposition if FFT supports transposed output*/
55 FFTW(plan) mpip[2];
56 #endif
57 MPI_Comm cart[2];
59 int N[3],M[3],K[3]; /*local length in transposed coordinate system (if not divisisable max)*/
60 int pN[3],pM[3], pK[3]; /*local length - not max but length for this processor*/
61 int oM[3],oK[3]; /*offset for current processor*/
62 int *iNin[3],*oNin[3],*iNout[3],*oNout[3]; /*size for each processor (if divisisable=max) for out(=split)
63 and in (=join) and offsets in transposed coordinate system*/
64 int C[3],rC[3]; /*global length (of the one global axes) */
65 /* C!=rC for real<->complex. then C=rC/2 but with potential padding*/
66 int P[2]; /*size of processor grid*/
67 /* int fftorder;*/
68 /* int direction;*/
69 /* int realcomplex;*/
70 int flags;
71 /*int N0,N1,M0,M1,K0,K1;*/
72 int NG,MG,KG;
73 /*int P[2];*/
74 int coor[2];
75 };
77 typedef struct fft5d_plan_t *fft5d_plan;
79 void fft5d_execute(fft5d_plan plan,fft5d_time times);
80 fft5d_plan fft5d_plan_3d(int N, int M, int K, MPI_Comm comm[2], int flags, t_complex** lin, t_complex** lin2);
81 void fft5d_local_size(fft5d_plan plan,int* N1,int* M0,int* K0,int* K1,int** coor);
82 void fft5d_destroy(fft5d_plan plan);
83 fft5d_plan fft5d_plan_3d_cart(int N, int M, int K, MPI_Comm comm, int P0, int flags, t_complex** lin, t_complex** lin2);
84 void fft5d_compare_data(const t_complex* lin, const t_complex* in, fft5d_plan plan, int bothLocal, int normarlize);
85 #endif /*FFTLIB_H_*/