SYCL: Avoid using no_init read accessor in rocFFT
[alexxy/gromacs.git] / src / gromacs / fft / fft5d.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2009-2017, The GROMACS development team.
5  * Copyright (c) 2019,2021, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36
37 #ifndef GMX_FFT_FFT5D_H
38 #define GMX_FFT_FFT5D_H
39
40 #include "config.h"
41
42 #ifdef NOGMX
43 /*#define GMX_MPI*/
44 /*#define GMX_FFT_FFTW3*/
45 FILE* debug;
46 #endif
47
48 #include "gromacs/fft/fft.h"
49 #include "gromacs/gpu_utils/hostallocator.h"
50 #include "gromacs/math/gmxcomplex.h"
51 #include "gromacs/utility/gmxmpi.h"
52
53 #if !GMX_MPI
54 double MPI_Wtime();
55 #endif
56
57 /*currently only special optimization for FFTE*/
58 #if GMX_FFT_FFTW3
59 #    include <fftw3.h>
60 #endif
61
62 #if !GMX_DOUBLE
63 #    define FFTW(x) fftwf_##x
64 #else
65 #    define FFTW(x) fftw_##x
66 #endif
67
68 #ifdef NOGMX
69 struct fft5d_time_t
70 {
71     double fft, local, mpi1, mpi2;
72 };
73 typedef struct fft5d_time_t* fft5d_time;
74 #else
75 #    include "gromacs/timing/wallcycle.h"
76 typedef gmx_wallcycle* fft5d_time;
77 #endif
78
79 namespace gmx
80 {
81 enum class PinningPolicy : int;
82 } // namespace gmx
83
84 typedef enum fft5d_flags_t
85 {
86     FFT5D_ORDER_YZ    = 1,
87     FFT5D_BACKWARD    = 2,
88     FFT5D_REALCOMPLEX = 4,
89     FFT5D_DEBUG       = 8,
90     FFT5D_NOMEASURE   = 16,
91     FFT5D_INPLACE     = 32,
92     FFT5D_NOMALLOC    = 64
93 } fft5d_flags;
94
95 struct fft5d_plan_t
96 {
97     t_complex* lin;
98     t_complex *lout, *lout2, *lout3;
99     gmx_fft_t* p1d[3]; /*1D plans*/
100 #if GMX_FFT_FFTW3
101     FFTW(plan) p2d; /*2D plan: used for 1D decomposition if FFT supports transposed output*/
102     FFTW(plan) p3d; /*3D plan: used for 0D decomposition if FFT supports transposed output*/
103     FFTW(plan) mpip[2];
104 #endif
105     MPI_Comm cart[2];
106
107     int  N[3], M[3], K[3]; /*local length in transposed coordinate system (if not divisisable max)*/
108     int  pN[3], pM[3], pK[3]; /*local length - not max but length for this processor*/
109     int  oM[3], oK[3];        /*offset for current processor*/
110     int *iNin[3], *oNin[3], *iNout[3],
111             *oNout[3]; /*size for each processor (if divisisable=max) for out(=split)
112                           and in (=join) and offsets in transposed coordinate system*/
113     int C[3], rC[3];   /*global length (of the one global axes) */
114     /* C!=rC for real<->complex. then C=rC/2 but with potential padding*/
115     int P[2]; /*size of processor grid*/
116               /*  int fftorder;*/
117               /*  int direction;*/
118               /*  int realcomplex;*/
119     int flags;
120     /*int N0,N1,M0,M1,K0,K1;*/
121     int NG, MG, KG;
122     /*int P[2];*/
123     int                coor[2];
124     int                nthreads;
125     gmx::PinningPolicy pinningPolicy;
126 };
127
128 typedef struct fft5d_plan_t* fft5d_plan;
129
130 void       fft5d_execute(fft5d_plan plan, int thread, fft5d_time times);
131 fft5d_plan fft5d_plan_3d(int         N,
132                          int         M,
133                          int         K,
134                          MPI_Comm    comm[2],
135                          int         flags,
136                          t_complex** lin,
137                          t_complex** lin2,
138                          t_complex** lout2,
139                          t_complex** lout3,
140                          int         nthreads,
141                          gmx::PinningPolicy realGridAllocationPinningPolicy = gmx::PinningPolicy::CannotBePinned);
142 void       fft5d_local_size(fft5d_plan plan, int* N1, int* M0, int* K0, int* K1, int** coor);
143 void       fft5d_destroy(fft5d_plan plan);
144 fft5d_plan fft5d_plan_3d_cart(int         N,
145                               int         M,
146                               int         K,
147                               MPI_Comm    comm,
148                               int         P0,
149                               int         flags,
150                               t_complex** lin,
151                               t_complex** lin2,
152                               t_complex** lout2,
153                               t_complex** lout3,
154                               int         nthreads);
155 void       fft5d_compare_data(const t_complex* lin, const t_complex* in, fft5d_plan plan, int bothLocal, int normarlize);
156
157 #endif