3349b63b79c189c1c6ac1d4a38d5466e904180a4
[alexxy/gromacs.git] / src / gromacs / ewald / pme-gpu-internal.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,2017, 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  *
38  * \brief This file contains internal function implementations
39  * for performing the PME calculations on GPU.
40  *
41  * \author Aleksei Iupinov <a.yupinov@gmail.com>
42  * \ingroup module_ewald
43  */
44
45 #include "gmxpre.h"
46
47 #include "pme-gpu-internal.h"
48
49 #include "config.h"
50
51 #include <list>
52 #include <string>
53
54 #include "gromacs/gpu_utils/gpu_utils.h"
55 #include "gromacs/math/invertmatrix.h"
56 #include "gromacs/math/units.h"
57 #include "gromacs/mdtypes/commrec.h"
58 #include "gromacs/utility/exceptions.h"
59 #include "gromacs/utility/gmxassert.h"
60 #include "gromacs/utility/logger.h"
61 #include "gromacs/utility/stringutil.h"
62
63 #include "pme-grid.h"
64 #include "pme-internal.h"
65
66 /*! \internal \brief
67  * Wrapper for getting a pointer to the plain C++ part of the GPU kernel parameters structure.
68  *
69  * \param[in] pmeGPU  The PME GPU structure.
70  * \returns The pointer to the kernel parameters.
71  */
72 static pme_gpu_kernel_params_base_t *pme_gpu_get_kernel_params_base_ptr(const pme_gpu_t *pmeGPU)
73 {
74     // reinterpret_cast is needed because the derived CUDA structure is not known in this file
75     auto *kernelParamsPtr = reinterpret_cast<pme_gpu_kernel_params_base_t *>(pmeGPU->kernelParams.get());
76     return kernelParamsPtr;
77 }
78
79 void pme_gpu_get_energy_virial(const pme_gpu_t *pmeGPU, real *energy, matrix virial)
80 {
81     GMX_ASSERT(energy, "Invalid energy output pointer in PME GPU");
82     unsigned int j = 0;
83     virial[XX][XX] = 0.25f * pmeGPU->staging.h_virialAndEnergy[j++];
84     virial[YY][YY] = 0.25f * pmeGPU->staging.h_virialAndEnergy[j++];
85     virial[ZZ][ZZ] = 0.25f * pmeGPU->staging.h_virialAndEnergy[j++];
86     virial[XX][YY] = virial[YY][XX] = 0.25f * pmeGPU->staging.h_virialAndEnergy[j++];
87     virial[XX][ZZ] = virial[ZZ][XX] = 0.25f * pmeGPU->staging.h_virialAndEnergy[j++];
88     virial[YY][ZZ] = virial[ZZ][YY] = 0.25f * pmeGPU->staging.h_virialAndEnergy[j++];
89     *energy        = 0.5f * pmeGPU->staging.h_virialAndEnergy[j++];
90 }
91
92 void pme_gpu_update_input_box(pme_gpu_t *pmeGPU, const matrix box)
93 {
94     auto        *kernelParamsPtr      = pme_gpu_get_kernel_params_base_ptr(pmeGPU);
95     kernelParamsPtr->step.boxVolume = box[XX][XX] * box[YY][YY] * box[ZZ][ZZ];
96     GMX_ASSERT(kernelParamsPtr->step.boxVolume != 0.0f, "Zero volume of the unit cell");
97
98 #if GMX_DOUBLE
99     GMX_THROW(gmx::NotImplementedError("PME is implemented for single-precision only on GPU"));
100 #else
101     matrix recipBox;
102     gmx::invertBoxMatrix(box, recipBox);
103     /* The GPU recipBox is transposed as compared to the CPU recipBox.
104      * Spread uses matrix columns (while solve and gather use rows).
105      * There is no particular reason for this; it might be further rethought/optimized for better access patterns.
106      */
107     const real newRecipBox[DIM][DIM] =
108     {
109         {recipBox[XX][XX], recipBox[YY][XX], recipBox[ZZ][XX]},
110         {             0.0, recipBox[YY][YY], recipBox[ZZ][YY]},
111         {             0.0,              0.0, recipBox[ZZ][ZZ]}
112     };
113     memcpy(kernelParamsPtr->step.recipBox, newRecipBox, sizeof(matrix));
114 #endif
115 }
116
117 void pme_gpu_start_step(pme_gpu_t *pmeGPU, bool needToUpdateBox, const matrix box, const rvec *h_coordinates)
118 {
119     pme_gpu_copy_input_coordinates(pmeGPU, h_coordinates);
120     if (needToUpdateBox)
121     {
122         pme_gpu_update_input_box(pmeGPU, box);
123     }
124 }
125
126 /*! \brief \libinternal
127  * The PME GPU reinitialization function that is called both at the end of any MD step and on any load balancing step.
128  *
129  * \param[in] pmeGPU            The PME GPU structure.
130  */
131 static void pme_gpu_reinit_step(const pme_gpu_t *pmeGPU)
132 {
133     pme_gpu_clear_grids(pmeGPU);
134     pme_gpu_clear_energy_virial(pmeGPU);
135 }
136
137 void pme_gpu_finish_step(const pme_gpu_t *pmeGPU, const bool bCalcF, const bool bCalcEnerVir)
138 {
139     /* Needed for copy back as well as timing events */
140     pme_gpu_synchronize(pmeGPU);
141
142     if (bCalcF && pme_gpu_performs_gather(pmeGPU))
143     {
144         pme_gpu_sync_output_forces(pmeGPU);
145     }
146     if (bCalcEnerVir && pme_gpu_performs_solve(pmeGPU))
147     {
148         pme_gpu_sync_output_energy_virial(pmeGPU);
149     }
150     pme_gpu_update_timings(pmeGPU);
151     pme_gpu_reinit_step(pmeGPU);
152 }
153
154 /*! \brief \libinternal
155  * (Re-)initializes all the PME GPU data related to the grid size and cut-off.
156  *
157  * \param[in] pmeGPU            The PME GPU structure.
158  */
159 static void pme_gpu_reinit_grids(pme_gpu_t *pmeGPU)
160 {
161     auto *kernelParamsPtr = pme_gpu_get_kernel_params_base_ptr(pmeGPU);
162     kernelParamsPtr->grid.ewaldFactor = (M_PI * M_PI) / (pmeGPU->common->ewaldcoeff_q * pmeGPU->common->ewaldcoeff_q);
163
164     /* The grid size variants */
165     for (int i = 0; i < DIM; i++)
166     {
167         kernelParamsPtr->grid.realGridSize[i]       = pmeGPU->common->nk[i];
168         kernelParamsPtr->grid.realGridSizeFP[i]     = (float)kernelParamsPtr->grid.realGridSize[i];
169         kernelParamsPtr->grid.realGridSizePadded[i] = kernelParamsPtr->grid.realGridSize[i];
170
171         // The complex grid currently uses no padding;
172         // if it starts to do so, then another test should be added for that
173         kernelParamsPtr->grid.complexGridSize[i]       = kernelParamsPtr->grid.realGridSize[i];
174         kernelParamsPtr->grid.complexGridSizePadded[i] = kernelParamsPtr->grid.realGridSize[i];
175     }
176     /* FFT: n real elements correspond to (n / 2 + 1) complex elements in minor dimension */
177     if (!pme_gpu_performs_FFT(pmeGPU))
178     {
179         // This allows for GPU spreading grid and CPU fftgrid to have the same layout, so that we can copy the data directly
180         kernelParamsPtr->grid.realGridSizePadded[ZZ] = (kernelParamsPtr->grid.realGridSize[ZZ] / 2 + 1) * 2;
181     }
182
183     /* GPU FFT: n real elements correspond to (n / 2 + 1) complex elements in minor dimension */
184     kernelParamsPtr->grid.complexGridSize[ZZ] /= 2;
185     kernelParamsPtr->grid.complexGridSize[ZZ]++;
186     kernelParamsPtr->grid.complexGridSizePadded[ZZ] = kernelParamsPtr->grid.complexGridSize[ZZ];
187
188     pme_gpu_realloc_and_copy_fract_shifts(pmeGPU);
189     pme_gpu_realloc_and_copy_bspline_values(pmeGPU);
190     pme_gpu_realloc_grids(pmeGPU);
191     pme_gpu_reinit_3dfft(pmeGPU);
192 }
193
194 /* Several GPU functions that refer to the CPU PME data live here.
195  * We would like to keep these away from the GPU-framework specific code for clarity,
196  * as well as compilation issues with MPI.
197  */
198
199 /*! \brief \libinternal
200  * Copies everything useful from the PME CPU to the PME GPU structure.
201  * The goal is to minimize interaction with the PME CPU structure in the GPU code.
202  *
203  * \param[in] pme         The PME structure.
204  */
205 static void pme_gpu_copy_common_data_from(const gmx_pme_t *pme)
206 {
207     /* TODO: Consider refactoring the CPU PME code to use the same structure,
208      * so that this function becomes 2 lines */
209     pme_gpu_t *pmeGPU             = pme->gpu;
210     pmeGPU->common->ngrids        = pme->ngrids;
211     pmeGPU->common->epsilon_r     = pme->epsilon_r;
212     pmeGPU->common->ewaldcoeff_q  = pme->ewaldcoeff_q;
213     pmeGPU->common->nk[XX]        = pme->nkx;
214     pmeGPU->common->nk[YY]        = pme->nky;
215     pmeGPU->common->nk[ZZ]        = pme->nkz;
216     pmeGPU->common->pmegrid_n[XX] = pme->pmegrid_nx;
217     pmeGPU->common->pmegrid_n[YY] = pme->pmegrid_ny;
218     pmeGPU->common->pmegrid_n[ZZ] = pme->pmegrid_nz;
219     pmeGPU->common->pme_order     = pme->pme_order;
220     for (int i = 0; i < DIM; i++)
221     {
222         pmeGPU->common->bsp_mod[i].assign(pme->bsp_mod[i], pme->bsp_mod[i] + pmeGPU->common->nk[i]);
223     }
224     const int cellCount = c_pmeNeighborUnitcellCount;
225     pmeGPU->common->fsh[XX].assign(pme->fshx, pme->fshx + cellCount * pme->nkx);
226     pmeGPU->common->fsh[YY].assign(pme->fshy, pme->fshy + cellCount * pme->nky);
227     pmeGPU->common->fsh[ZZ].assign(pme->fshz, pme->fshz + cellCount * pme->nkz);
228     pmeGPU->common->nn[XX].assign(pme->nnx, pme->nnx + cellCount * pme->nkx);
229     pmeGPU->common->nn[YY].assign(pme->nny, pme->nny + cellCount * pme->nky);
230     pmeGPU->common->nn[ZZ].assign(pme->nnz, pme->nnz + cellCount * pme->nkz);
231     pmeGPU->common->runMode = pme->runMode;
232 }
233
234 /*! \brief \libinternal
235  * Finds out if PME with given inputs is possible to run on GPU.
236  *
237  * \param[in]  pme          The PME structure.
238  * \param[out] error        The error message if the input is not supported on GPU.
239  * \returns                 True if this PME input is possible to run on GPU, false otherwise.
240  */
241 static bool pme_gpu_check_restrictions(const gmx_pme_t *pme, std::string *error)
242 {
243     std::list<std::string> errorReasons;
244     if (pme->nnodes != 1)
245     {
246         errorReasons.push_back("PME decomposition");
247     }
248     if (pme->pme_order != 4)
249     {
250         errorReasons.push_back("interpolation orders other than 4");
251     }
252     if (pme->bFEP)
253     {
254         errorReasons.push_back("free energy calculations (multiple grids)");
255     }
256     if (pme->doLJ)
257     {
258         errorReasons.push_back("Lennard-Jones PME");
259     }
260 #if GMX_DOUBLE
261     {
262         errorReasons.push_back("double precision");
263     }
264 #endif
265 #if GMX_GPU != GMX_GPU_CUDA
266     {
267         errorReasons.push_back("non-CUDA build of GROMACS");
268     }
269 #endif
270
271     bool inputSupported = errorReasons.empty();
272     if (!inputSupported && error)
273     {
274         std::string regressionTestMarker = "PME GPU does not support";
275         // this prefix is tested for in the regression tests script gmxtest.pl
276         *error = regressionTestMarker + ": " + gmx::joinStrings(errorReasons, "; ") + ".";
277     }
278     return inputSupported;
279 }
280
281 /*! \libinternal \brief
282  * Initializes the PME GPU data at the beginning of the run.
283  *
284  * \param[in,out] pme       The PME structure.
285  * \param[in,out] gpuInfo   The GPU information structure.
286  * \param[in]     mdlog     The logger.
287  * \param[in]     cr        The communication structure.
288  */
289 static void pme_gpu_init(gmx_pme_t *pme, gmx_device_info_t *gpuInfo, const gmx::MDLogger &mdlog,
290                          const t_commrec *cr)
291 {
292     std::string errorString;
293     bool        canRunOnGpu = pme_gpu_check_restrictions(pme, &errorString);
294     if (!canRunOnGpu)
295     {
296         GMX_THROW(gmx::NotImplementedError(errorString));
297     }
298
299     pme->gpu          = new pme_gpu_t();
300     pme_gpu_t *pmeGPU = pme->gpu;
301     pmeGPU->common = std::shared_ptr<pme_shared_t>(new pme_shared_t());
302
303     /* These settings are set here for the whole run; dynamic ones are set in pme_gpu_reinit() */
304     /* A convenience variable. */
305     pmeGPU->settings.useDecomposition = (pme->nnodes == 1);
306     /* TODO: CPU gather with GPU spread is broken due to different theta/dtheta layout. */
307     pmeGPU->settings.performGPUGather = true;
308
309     // GPU initialization
310     init_gpu(mdlog, cr->nodeid, gpuInfo);
311     pmeGPU->deviceInfo = gpuInfo;
312
313     pme_gpu_init_internal(pmeGPU);
314     pme_gpu_init_sync_events(pmeGPU);
315     pme_gpu_alloc_energy_virial(pmeGPU);
316
317     pme_gpu_copy_common_data_from(pme);
318
319     GMX_ASSERT(pmeGPU->common->epsilon_r != 0.0f, "PME GPU: bad electrostatic coefficient");
320
321     auto *kernelParamsPtr = pme_gpu_get_kernel_params_base_ptr(pmeGPU);
322     kernelParamsPtr->constants.elFactor = ONE_4PI_EPS0 / pmeGPU->common->epsilon_r;
323 }
324
325 void pme_gpu_reinit(gmx_pme_t *pme, gmx_device_info_t *gpuInfo, const gmx::MDLogger &mdlog, const t_commrec *cr)
326 {
327     if (!pme_gpu_active(pme))
328     {
329         return;
330     }
331
332     if (!pme->gpu)
333     {
334         /* First-time initialization */
335         pme_gpu_init(pme, gpuInfo, mdlog, cr);
336     }
337     else
338     {
339         /* After this call nothing in the GPU code should refer to the gmx_pme_t *pme itself - until the next pme_gpu_reinit */
340         pme_gpu_copy_common_data_from(pme);
341     }
342     /* GPU FFT will only get used for a single rank.*/
343     pme->gpu->settings.performGPUFFT   = (pme->gpu->common->runMode == PmeRunMode::GPU) && !pme_gpu_uses_dd(pme->gpu);
344     pme->gpu->settings.performGPUSolve = (pme->gpu->common->runMode == PmeRunMode::GPU);
345
346     pme_gpu_reinit_grids(pme->gpu);
347     pme_gpu_reinit_step(pme->gpu);
348 }
349
350 void pme_gpu_destroy(pme_gpu_t *pmeGPU)
351 {
352     /* Free lots of data */
353     pme_gpu_free_energy_virial(pmeGPU);
354     pme_gpu_free_bspline_values(pmeGPU);
355     pme_gpu_free_forces(pmeGPU);
356     pme_gpu_free_coordinates(pmeGPU);
357     pme_gpu_free_coefficients(pmeGPU);
358     pme_gpu_free_spline_data(pmeGPU);
359     pme_gpu_free_grid_indices(pmeGPU);
360     pme_gpu_free_fract_shifts(pmeGPU);
361     pme_gpu_free_grids(pmeGPU);
362
363     pme_gpu_destroy_3dfft(pmeGPU);
364     pme_gpu_destroy_sync_events(pmeGPU);
365
366     /* Free the GPU-framework specific data last */
367     pme_gpu_destroy_specific(pmeGPU);
368
369     delete pmeGPU;
370 }
371
372 void pme_gpu_reinit_atoms(pme_gpu_t *pmeGPU, const int nAtoms, const real *charges)
373 {
374     auto      *kernelParamsPtr = pme_gpu_get_kernel_params_base_ptr(pmeGPU);
375     kernelParamsPtr->atoms.nAtoms = nAtoms;
376     const int  alignment = pme_gpu_get_atom_data_alignment(pmeGPU);
377     pmeGPU->nAtomsPadded = ((nAtoms + alignment - 1) / alignment) * alignment;
378     const int  nAtomsAlloc   = c_usePadding ? pmeGPU->nAtomsPadded : nAtoms;
379     const bool haveToRealloc = (pmeGPU->nAtomsAlloc < nAtomsAlloc); /* This check might be redundant, but is logical */
380     pmeGPU->nAtomsAlloc = nAtomsAlloc;
381
382 #if GMX_DOUBLE
383     GMX_RELEASE_ASSERT(false, "Only single precision supported");
384     GMX_UNUSED_VALUE(charges);
385 #else
386     pme_gpu_realloc_and_copy_input_coefficients(pmeGPU, reinterpret_cast<const float *>(charges));
387     /* Could also be checked for haveToRealloc, but the copy always needs to be performed */
388 #endif
389
390     if (haveToRealloc)
391     {
392         pme_gpu_realloc_coordinates(pmeGPU);
393         pme_gpu_realloc_forces(pmeGPU);
394         pme_gpu_realloc_spline_data(pmeGPU);
395         pme_gpu_realloc_grid_indices(pmeGPU);
396     }
397 }