Test gmxapi and clients through Py 3.9.
[alexxy/gromacs.git] / src / gromacs / fft / gpu_3dfft.cu
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,2017,2018,2019,2020,2021, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35
36 /*! \internal \file
37  *  \brief Implements GPU 3D FFT routines for CUDA.
38  *
39  *  \author Aleksei Iupinov <a.yupinov@gmail.com>
40  *  \author Mark Abraham <mark.j.abraham@gmail.com>
41  *  \ingroup module_fft
42  */
43
44 #include "gmxpre.h"
45
46 #include "gpu_3dfft.h"
47
48 #include <cufft.h>
49
50 #include "gromacs/gpu_utils/device_stream.h"
51 #include "gromacs/utility/fatalerror.h"
52 #include "gromacs/utility/gmxassert.h"
53
54 namespace gmx
55 {
56
57 class Gpu3dFft::Impl
58 {
59 public:
60     Impl(ivec                 realGridSize,
61          ivec                 realGridSizePadded,
62          ivec                 complexGridSizePadded,
63          bool                 useDecomposition,
64          bool                 performOutOfPlaceFFT,
65          const DeviceContext& context,
66          const DeviceStream&  pmeStream,
67          DeviceBuffer<float>  realGrid,
68          DeviceBuffer<float>  complexGrid);
69     ~Impl();
70
71     cufftHandle   planR2C_;
72     cufftHandle   planC2R_;
73     cufftReal*    realGrid_;
74     cufftComplex* complexGrid_;
75 };
76
77 static void handleCufftError(cufftResult_t status, const char* msg)
78 {
79     if (status != CUFFT_SUCCESS)
80     {
81         gmx_fatal(FARGS, "%s (error code %d)\n", msg, status);
82     }
83 }
84
85 Gpu3dFft::Impl::Impl(ivec       realGridSize,
86                      ivec       realGridSizePadded,
87                      ivec       complexGridSizePadded,
88                      const bool useDecomposition,
89                      const bool /*performOutOfPlaceFFT*/,
90                      const DeviceContext& /*context*/,
91                      const DeviceStream& pmeStream,
92                      DeviceBuffer<float> realGrid,
93                      DeviceBuffer<float> complexGrid) :
94     realGrid_(reinterpret_cast<cufftReal*>(realGrid)),
95     complexGrid_(reinterpret_cast<cufftComplex*>(complexGrid))
96 {
97     GMX_RELEASE_ASSERT(!useDecomposition, "FFT decomposition not implemented");
98
99     const int complexGridSizePaddedTotal =
100             complexGridSizePadded[XX] * complexGridSizePadded[YY] * complexGridSizePadded[ZZ];
101     const int realGridSizePaddedTotal =
102             realGridSizePadded[XX] * realGridSizePadded[YY] * realGridSizePadded[ZZ];
103
104     GMX_RELEASE_ASSERT(realGrid_, "Bad (null) input real-space grid");
105     GMX_RELEASE_ASSERT(complexGrid_, "Bad (null) input complex grid");
106
107     cufftResult_t result;
108     /* Commented code for a simple 3D grid with no padding */
109     /*
110        result = cufftPlan3d(&planR2C_, realGridSize[XX], realGridSize[YY], realGridSize[ZZ],
111        CUFFT_R2C); handleCufftError(result, "cufftPlan3d R2C plan failure");
112
113        result = cufftPlan3d(&planC2R_, realGridSize[XX], realGridSize[YY], realGridSize[ZZ],
114        CUFFT_C2R); handleCufftError(result, "cufftPlan3d C2R plan failure");
115      */
116
117     const int rank = 3, batch = 1;
118     result = cufftPlanMany(&planR2C_,
119                            rank,
120                            realGridSize,
121                            realGridSizePadded,
122                            1,
123                            realGridSizePaddedTotal,
124                            complexGridSizePadded,
125                            1,
126                            complexGridSizePaddedTotal,
127                            CUFFT_R2C,
128                            batch);
129     handleCufftError(result, "cufftPlanMany R2C plan failure");
130
131     result = cufftPlanMany(&planC2R_,
132                            rank,
133                            realGridSize,
134                            complexGridSizePadded,
135                            1,
136                            complexGridSizePaddedTotal,
137                            realGridSizePadded,
138                            1,
139                            realGridSizePaddedTotal,
140                            CUFFT_C2R,
141                            batch);
142     handleCufftError(result, "cufftPlanMany C2R plan failure");
143
144     cudaStream_t stream = pmeStream.stream();
145     GMX_RELEASE_ASSERT(stream, "Can not use the default CUDA stream for PME cuFFT");
146
147     result = cufftSetStream(planR2C_, stream);
148     handleCufftError(result, "cufftSetStream R2C failure");
149
150     result = cufftSetStream(planC2R_, stream);
151     handleCufftError(result, "cufftSetStream C2R failure");
152 }
153
154 Gpu3dFft::Impl::~Impl()
155 {
156     cufftResult_t result;
157     result = cufftDestroy(planR2C_);
158     handleCufftError(result, "cufftDestroy R2C failure");
159     result = cufftDestroy(planC2R_);
160     handleCufftError(result, "cufftDestroy C2R failure");
161 }
162
163 void Gpu3dFft::perform3dFft(gmx_fft_direction dir, CommandEvent* /*timingEvent*/)
164 {
165     cufftResult_t result;
166     if (dir == GMX_FFT_REAL_TO_COMPLEX)
167     {
168         result = cufftExecR2C(impl_->planR2C_, impl_->realGrid_, impl_->complexGrid_);
169         handleCufftError(result, "cuFFT R2C execution failure");
170     }
171     else
172     {
173         result = cufftExecC2R(impl_->planC2R_, impl_->complexGrid_, impl_->realGrid_);
174         handleCufftError(result, "cuFFT C2R execution failure");
175     }
176 }
177
178 Gpu3dFft::Gpu3dFft(ivec                 realGridSize,
179                    ivec                 realGridSizePadded,
180                    ivec                 complexGridSizePadded,
181                    const bool           useDecomposition,
182                    const bool           performOutOfPlaceFFT,
183                    const DeviceContext& context,
184                    const DeviceStream&  pmeStream,
185                    DeviceBuffer<float>  realGrid,
186                    DeviceBuffer<float>  complexGrid) :
187     impl_(std::make_unique<Impl>(realGridSize,
188                                  realGridSizePadded,
189                                  complexGridSizePadded,
190                                  useDecomposition,
191                                  performOutOfPlaceFFT,
192                                  context,
193                                  pmeStream,
194                                  realGrid,
195                                  complexGrid))
196 {
197 }
198
199 Gpu3dFft::~Gpu3dFft() = default;
200
201 } // namespace gmx