e22fe645928d5555b43087003e71c987920dd866
[alexxy/gromacs.git] / src / gromacs / ewald / pme_gpu_internal.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,2017,2018,2019, 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 definitions for performing the PME calculations on GPU.
39  * These are not meant to be exposed outside of the PME GPU code.
40  * As of now, their bodies are still in the common pme_gpu.cpp files.
41  *
42  * \author Aleksei Iupinov <a.yupinov@gmail.com>
43  * \ingroup module_ewald
44  */
45
46 #ifndef GMX_EWALD_PME_GPU_INTERNAL_H
47 #define GMX_EWALD_PME_GPU_INTERNAL_H
48
49 #include "gromacs/fft/fft.h"                   // for the gmx_fft_direction enum
50 #include "gromacs/gpu_utils/gpu_macros.h"      // for the GPU_FUNC_ macros
51 #include "gromacs/utility/arrayref.h"
52
53 #include "pme_gpu_types_host.h"                     // for the inline functions accessing PmeGpu members
54
55 struct gmx_hw_info_t;
56 struct gmx_gpu_opt_t;
57 struct gmx_pme_t;                              // only used in pme_gpu_reinit
58 struct gmx_wallclock_gpu_pme_t;
59 class PmeAtomComm;
60 struct t_complex;
61
62 namespace gmx
63 {
64 class MDLogger;
65 }
66
67 //! Type of spline data
68 enum class PmeSplineDataType
69 {
70     Values,      // theta
71     Derivatives, // dtheta
72 };               //TODO move this into new and shiny pme.h (pme-types.h?)
73
74 //! PME grid dimension ordering (from major to minor)
75 enum class GridOrdering
76 {
77     YZX,
78     XYZ
79 };
80
81 /*! \libinternal \brief
82  * Returns the number of atoms per chunk in the atom charges/coordinates data layout.
83  * Depends on CUDA-specific block sizes, needed for the atom data padding.
84  *
85  * \param[in] pmeGpu            The PME GPU structure.
86  * \returns   Number of atoms in a single GPU atom data chunk.
87  */
88 int pme_gpu_get_atom_data_alignment(const PmeGpu *pmeGpu);
89
90 /*! \libinternal \brief
91  * Returns the number of atoms per chunk in the atom spline theta/dtheta data layout.
92  *
93  * \param[in] pmeGpu            The PME GPU structure.
94  * \returns   Number of atoms in a single GPU atom spline data chunk.
95  */
96 int pme_gpu_get_atoms_per_warp(const PmeGpu *pmeGpu);
97
98 /*! \libinternal \brief
99  * Synchronizes the current computation, waiting for the GPU kernels/transfers to finish.
100  *
101  * \param[in] pmeGpu            The PME GPU structure.
102  */
103 GPU_FUNC_QUALIFIER void pme_gpu_synchronize(const PmeGpu *GPU_FUNC_ARGUMENT(pmeGpu)) GPU_FUNC_TERM;
104
105 /*! \libinternal \brief
106  * Allocates the fixed size energy and virial buffer both on GPU and CPU.
107  *
108  * \param[in,out] pmeGpu            The PME GPU structure.
109  */
110 void pme_gpu_alloc_energy_virial(PmeGpu *pmeGpu);
111
112 /*! \libinternal \brief
113  * Frees the energy and virial memory both on GPU and CPU.
114  *
115  * \param[in] pmeGpu            The PME GPU structure.
116  */
117 void pme_gpu_free_energy_virial(PmeGpu *pmeGpu);
118
119 /*! \libinternal \brief
120  * Clears the energy and virial memory on GPU with 0.
121  * Should be called at the end of PME computation which returned energy/virial.
122  *
123  * \param[in] pmeGpu            The PME GPU structure.
124  */
125 void pme_gpu_clear_energy_virial(const PmeGpu *pmeGpu);
126
127 /*! \libinternal \brief
128  * Reallocates and copies the pre-computed B-spline values to the GPU.
129  *
130  * \param[in,out] pmeGpu             The PME GPU structure.
131  */
132 void pme_gpu_realloc_and_copy_bspline_values(PmeGpu *pmeGpu);
133
134 /*! \libinternal \brief
135  * Frees the pre-computed B-spline values on the GPU (and the transfer CPU buffers).
136  *
137  * \param[in] pmeGpu             The PME GPU structure.
138  */
139 void pme_gpu_free_bspline_values(const PmeGpu *pmeGpu);
140
141 /*! \libinternal \brief
142  * Reallocates the GPU buffer for the PME forces.
143  *
144  * \param[in] pmeGpu             The PME GPU structure.
145  */
146 void pme_gpu_realloc_forces(PmeGpu *pmeGpu);
147
148 /*! \libinternal \brief
149  * Frees the GPU buffer for the PME forces.
150  *
151  * \param[in] pmeGpu             The PME GPU structure.
152  */
153 void pme_gpu_free_forces(const PmeGpu *pmeGpu);
154
155 /*! \libinternal \brief
156  * Copies the forces from the CPU buffer to the GPU (to reduce them with the PME GPU gathered forces).
157  * To be called e.g. after the bonded calculations.
158  *
159  * \param[in] pmeGpu             The PME GPU structure.
160  */
161 void pme_gpu_copy_input_forces(PmeGpu *pmeGpu);
162
163 /*! \libinternal \brief
164  * Copies the forces from the GPU to the CPU buffer. To be called after the gathering stage.
165  *
166  * \param[in] pmeGpu             The PME GPU structure.
167  */
168 void pme_gpu_copy_output_forces(PmeGpu *pmeGpu);
169
170 /*! \libinternal \brief
171  * Checks whether work in the PME GPU stream has completed.
172  *
173  * \param[in] pmeGpu            The PME GPU structure.
174  *
175  * \returns                     True if work in the PME stream has completed.
176  */
177 bool pme_gpu_stream_query(const PmeGpu *pmeGpu);
178
179 /*! \libinternal \brief
180  * Reallocates the input coordinates buffer on the GPU (and clears the padded part if needed).
181  *
182  * \param[in] pmeGpu            The PME GPU structure.
183  *
184  * Needs to be called on every DD step/in the beginning.
185  */
186 void pme_gpu_realloc_coordinates(const PmeGpu *pmeGpu);
187
188 /*! \libinternal \brief
189  * Copies the input coordinates from the CPU buffer onto the GPU.
190  *
191  * \param[in] pmeGpu            The PME GPU structure.
192  * \param[in] h_coordinates     Input coordinates (XYZ rvec array).
193  *
194  * Needs to be called for every PME computation. The coordinates are then used in the spline calculation.
195  */
196 GPU_FUNC_QUALIFIER void pme_gpu_copy_input_coordinates(const PmeGpu *GPU_FUNC_ARGUMENT(pmeGpu),
197                                                        const rvec   *GPU_FUNC_ARGUMENT(h_coordinates)) GPU_FUNC_TERM;
198
199 /*! \libinternal \brief
200  * Frees the coordinates on the GPU.
201  *
202  * \param[in] pmeGpu            The PME GPU structure.
203  */
204 void pme_gpu_free_coordinates(const PmeGpu *pmeGpu);
205
206 /*! \libinternal \brief
207  * Reallocates the buffer on the GPU and copies the charges/coefficients from the CPU buffer.
208  * Clears the padded part if needed.
209  *
210  * \param[in] pmeGpu            The PME GPU structure.
211  * \param[in] h_coefficients    The input atom charges/coefficients.
212  *
213  * Does not need to be done for every PME computation, only whenever the local charges change.
214  * (So, in the beginning of the run, or on DD step).
215  */
216 void pme_gpu_realloc_and_copy_input_coefficients(const PmeGpu    *pmeGpu,
217                                                  const float     *h_coefficients);
218
219 /*! \libinternal \brief
220  * Frees the charges/coefficients on the GPU.
221  *
222  * \param[in] pmeGpu             The PME GPU structure.
223  */
224 void pme_gpu_free_coefficients(const PmeGpu *pmeGpu);
225
226 /*! \libinternal \brief
227  * Reallocates the buffers on the GPU and the host for the atoms spline data.
228  *
229  * \param[in,out] pmeGpu            The PME GPU structure.
230  */
231 void pme_gpu_realloc_spline_data(PmeGpu *pmeGpu);
232
233 /*! \libinternal \brief
234  * Frees the buffers on the GPU for the atoms spline data.
235  *
236  * \param[in] pmeGpu            The PME GPU structure.
237  */
238 void pme_gpu_free_spline_data(const PmeGpu *pmeGpu);
239
240 /*! \libinternal \brief
241  * Reallocates the buffers on the GPU and the host for the particle gridline indices.
242  *
243  * \param[in,out] pmeGpu            The PME GPU structure.
244  */
245 void pme_gpu_realloc_grid_indices(PmeGpu *pmeGpu);
246
247 /*! \libinternal \brief
248  * Frees the buffer on the GPU for the particle gridline indices.
249  *
250  * \param[in] pmeGpu            The PME GPU structure.
251  */
252 void pme_gpu_free_grid_indices(const PmeGpu *pmeGpu);
253
254 /*! \libinternal \brief
255  * Reallocates the real space grid and the complex reciprocal grid (if needed) on the GPU.
256  *
257  * \param[in] pmeGpu            The PME GPU structure.
258  */
259 void pme_gpu_realloc_grids(PmeGpu *pmeGpu);
260
261 /*! \libinternal \brief
262  * Frees the real space grid and the complex reciprocal grid (if needed) on the GPU.
263  *
264  * \param[in] pmeGpu            The PME GPU structure.
265  */
266 void pme_gpu_free_grids(const PmeGpu *pmeGpu);
267
268 /*! \libinternal \brief
269  * Clears the real space grid on the GPU.
270  * Should be called at the end of each computation.
271  *
272  * \param[in] pmeGpu            The PME GPU structure.
273  */
274 void pme_gpu_clear_grids(const PmeGpu *pmeGpu);
275
276 /*! \libinternal \brief
277  * Reallocates and copies the pre-computed fractional coordinates' shifts to the GPU.
278  *
279  * \param[in] pmeGpu            The PME GPU structure.
280  */
281 void pme_gpu_realloc_and_copy_fract_shifts(PmeGpu *pmeGpu);
282
283 /*! \libinternal \brief
284  * Frees the pre-computed fractional coordinates' shifts on the GPU.
285  *
286  * \param[in] pmeGpu            The PME GPU structure.
287  */
288 void pme_gpu_free_fract_shifts(const PmeGpu *pmeGpu);
289
290 /*! \libinternal \brief
291  * Copies the input real-space grid from the host to the GPU.
292  *
293  * \param[in] pmeGpu   The PME GPU structure.
294  * \param[in] h_grid   The host-side grid buffer.
295  */
296 void pme_gpu_copy_input_gather_grid(const PmeGpu *pmeGpu,
297                                     float        *h_grid);
298
299 /*! \libinternal \brief
300  * Copies the output real-space grid from the GPU to the host.
301  *
302  * \param[in] pmeGpu   The PME GPU structure.
303  * \param[out] h_grid  The host-side grid buffer.
304  */
305 void pme_gpu_copy_output_spread_grid(const PmeGpu *pmeGpu,
306                                      float        *h_grid);
307
308 /*! \libinternal \brief
309  * Copies the spread output spline data and gridline indices from the GPU to the host.
310  *
311  * \param[in] pmeGpu   The PME GPU structure.
312  */
313 void pme_gpu_copy_output_spread_atom_data(const PmeGpu *pmeGpu);
314
315 /*! \libinternal \brief
316  * Copies the gather input spline data and gridline indices from the host to the GPU.
317  *
318  * \param[in] pmeGpu   The PME GPU structure.
319  */
320 void pme_gpu_copy_input_gather_atom_data(const PmeGpu *pmeGpu);
321
322 /*! \libinternal \brief
323  * Waits for the grid copying to the host-side buffer after spreading to finish.
324  *
325  * \param[in] pmeGpu  The PME GPU structure.
326  */
327 void pme_gpu_sync_spread_grid(const PmeGpu *pmeGpu);
328
329 /*! \libinternal \brief
330  * Does the one-time GPU-framework specific PME initialization.
331  * For CUDA, the PME stream is created with the highest priority.
332  *
333  * \param[in] pmeGpu  The PME GPU structure.
334  */
335 void pme_gpu_init_internal(PmeGpu *pmeGpu);
336
337 /*! \libinternal \brief
338  * Destroys the PME GPU-framework specific data.
339  * Should be called last in the PME GPU destructor.
340  *
341  * \param[in] pmeGpu  The PME GPU structure.
342  */
343 void pme_gpu_destroy_specific(const PmeGpu *pmeGpu);
344
345 /*! \libinternal \brief
346  * Initializes the CUDA FFT structures.
347  *
348  * \param[in] pmeGpu  The PME GPU structure.
349  */
350 void pme_gpu_reinit_3dfft(const PmeGpu *pmeGpu);
351
352 /*! \libinternal \brief
353  * Destroys the CUDA FFT structures.
354  *
355  * \param[in] pmeGpu  The PME GPU structure.
356  */
357 void pme_gpu_destroy_3dfft(const PmeGpu *pmeGpu);
358
359 /* Several GPU event-based timing functions that live in pme_gpu_timings.cpp */
360
361 /*! \libinternal \brief
362  * Finalizes all the active PME GPU stage timings for the current computation. Should be called at the end of every computation.
363  *
364  * \param[in] pmeGpu         The PME GPU structure.
365  */
366 void pme_gpu_update_timings(const PmeGpu *pmeGpu);
367
368 /*! \libinternal \brief
369  * Updates the internal list of active PME GPU stages (if timings are enabled).
370  *
371  * \param[in] pmeGpu         The PME GPU data structure.
372  */
373 void pme_gpu_reinit_timings(const PmeGpu *pmeGpu);
374
375 /*! \brief
376  * Resets the PME GPU timings. To be called at the reset MD step.
377  *
378  * \param[in] pmeGpu         The PME GPU structure.
379  */
380 void pme_gpu_reset_timings(const PmeGpu *pmeGpu);
381
382 /*! \libinternal \brief
383  * Copies the PME GPU timings to the gmx_wallclock_gpu_t structure (for log output). To be called at the run end.
384  *
385  * \param[in] pmeGpu         The PME GPU structure.
386  * \param[in] timings        The gmx_wallclock_gpu_pme_t structure.
387  */
388 void pme_gpu_get_timings(const PmeGpu            *pmeGpu,
389                          gmx_wallclock_gpu_pme_t *timings);
390
391 /* The PME stages themselves */
392
393 /*! \libinternal \brief
394  * A GPU spline computation and charge spreading function.
395  *
396  * \param[in]  pmeGpu          The PME GPU structure.
397  * \param[in]  gridIndex       Index of the PME grid - unused, assumed to be 0.
398  * \param[out] h_grid          The host-side grid buffer (used only if the result of the spread is expected on the host,
399  *                             e.g. testing or host-side FFT)
400  * \param[in]  computeSplines  Should the computation of spline parameters and gridline indices be performed.
401  * \param[in]  spreadCharges   Should the charges/coefficients be spread on the grid.
402  */
403 GPU_FUNC_QUALIFIER void pme_gpu_spread(const PmeGpu    *GPU_FUNC_ARGUMENT(pmeGpu),
404                                        int              GPU_FUNC_ARGUMENT(gridIndex),
405                                        real            *GPU_FUNC_ARGUMENT(h_grid),
406                                        bool             GPU_FUNC_ARGUMENT(computeSplines),
407                                        bool             GPU_FUNC_ARGUMENT(spreadCharges)) GPU_FUNC_TERM;
408
409 /*! \libinternal \brief
410  * 3D FFT R2C/C2R routine.
411  *
412  * \param[in]  pmeGpu          The PME GPU structure.
413  * \param[in]  direction       Transform direction (real-to-complex or complex-to-real)
414  * \param[in]  gridIndex       Index of the PME grid - unused, assumed to be 0.
415  */
416 void pme_gpu_3dfft(const PmeGpu          *pmeGpu,
417                    enum gmx_fft_direction direction,
418                    int                    gridIndex);
419
420 /*! \libinternal \brief
421  * A GPU Fourier space solving function.
422  *
423  * \param[in]     pmeGpu                  The PME GPU structure.
424  * \param[in,out] h_grid                  The host-side input and output Fourier grid buffer (used only with testing or host-side FFT)
425  * \param[in]     gridOrdering            Specifies the dimenion ordering of the complex grid. TODO: store this information?
426  * \param[in]     computeEnergyAndVirial  Tells if the energy and virial computation should also be performed.
427  */
428 GPU_FUNC_QUALIFIER void pme_gpu_solve(const PmeGpu    *GPU_FUNC_ARGUMENT(pmeGpu),
429                                       t_complex       *GPU_FUNC_ARGUMENT(h_grid),
430                                       GridOrdering     GPU_FUNC_ARGUMENT(gridOrdering),
431                                       bool             GPU_FUNC_ARGUMENT(computeEnergyAndVirial)) GPU_FUNC_TERM;
432
433 /*! \libinternal \brief
434  * A GPU force gathering function.
435  *
436  * \param[in]     pmeGpu           The PME GPU structure.
437  * \param[in]     forceTreatment   Tells how data in h_forces should be treated.
438  *                                 TODO: determine efficiency/balance of host/device-side reductions.
439  * \param[in]     h_grid           The host-side grid buffer (used only in testing mode)
440  */
441 GPU_FUNC_QUALIFIER void pme_gpu_gather(PmeGpu                *GPU_FUNC_ARGUMENT(pmeGpu),
442                                        PmeForceOutputHandling GPU_FUNC_ARGUMENT(forceTreatment),
443                                        const float           *GPU_FUNC_ARGUMENT(h_grid)) GPU_FUNC_TERM;
444
445 /*! \brief Return pointer to device copy of coordinate data. */
446 GPU_FUNC_QUALIFIER void * pme_gpu_get_kernelparam_coordinates(const PmeGpu *GPU_FUNC_ARGUMENT(pmeGpu)) GPU_FUNC_TERM_WITH_RETURN(nullptr);
447
448 /* The inlined convenience PME GPU status getters */
449
450 /*! \libinternal \brief
451  * Tells if PME runs on multiple GPUs with the decomposition.
452  *
453  * \param[in] pmeGpu         The PME GPU structure.
454  * \returns                  True if PME runs on multiple GPUs, false otherwise.
455  */
456 inline bool pme_gpu_uses_dd(const PmeGpu *pmeGpu)
457 {
458     return !pmeGpu->settings.useDecomposition;
459 }
460
461 /*! \libinternal \brief
462  * Tells if PME performs the gathering stage on GPU.
463  *
464  * \param[in] pmeGpu         The PME GPU structure.
465  * \returns                  True if the gathering is performed on GPU, false otherwise.
466  */
467 inline bool pme_gpu_performs_gather(const PmeGpu *pmeGpu)
468 {
469     return pmeGpu->settings.performGPUGather;
470 }
471
472 /*! \libinternal \brief
473  * Tells if PME performs the FFT stages on GPU.
474  *
475  * \param[in] pmeGpu         The PME GPU structure.
476  * \returns                  True if FFT is performed on GPU, false otherwise.
477  */
478 inline bool pme_gpu_performs_FFT(const PmeGpu *pmeGpu)
479 {
480     return pmeGpu->settings.performGPUFFT;
481 }
482
483 /*! \libinternal \brief
484  * Tells if PME performs the grid (un-)wrapping on GPU.
485  *
486  * \param[in] pmeGpu         The PME GPU structure.
487  * \returns                  True if (un-)wrapping is performed on GPU, false otherwise.
488  */
489 inline bool pme_gpu_performs_wrapping(const PmeGpu *pmeGpu)
490 {
491     return pmeGpu->settings.useDecomposition;
492 }
493
494 /*! \libinternal \brief
495  * Tells if PME performs the grid solving on GPU.
496  *
497  * \param[in] pmeGpu         The PME GPU structure.
498  * \returns                  True if solving is performed on GPU, false otherwise.
499  */
500 inline bool pme_gpu_performs_solve(const PmeGpu *pmeGpu)
501 {
502     return pmeGpu->settings.performGPUSolve;
503 }
504
505 /*! \libinternal \brief
506  * Enables or disables the testing mode.
507  * Testing mode only implies copying all the outputs, even the intermediate ones, to the host,
508  * and also makes the copies synchronous.
509  *
510  * \param[in] pmeGpu             The PME GPU structure.
511  * \param[in] testing            Should the testing mode be enabled, or disabled.
512  */
513 inline void pme_gpu_set_testing(PmeGpu *pmeGpu, bool testing)
514 {
515     if (pmeGpu)
516     {
517         pmeGpu->settings.copyAllOutputs = testing;
518         pmeGpu->settings.transferKind   = testing ? GpuApiCallBehavior::Sync : GpuApiCallBehavior::Async;
519     }
520 }
521
522 /*! \libinternal \brief
523  * Tells if PME is in the testing mode.
524  *
525  * \param[in] pmeGpu             The PME GPU structure.
526  * \returns                      true if testing mode is enabled, false otherwise.
527  */
528 inline bool pme_gpu_is_testing(const PmeGpu *pmeGpu)
529 {
530     return pmeGpu->settings.copyAllOutputs;
531 }
532
533 /* A block of C++ functions that live in pme_gpu_internal.cpp */
534
535 /*! \libinternal \brief
536  * Returns the energy and virial GPU outputs, useful for testing.
537  *
538  * It is the caller's responsibility to be aware of whether the GPU
539  * handled the solve stage.
540  *
541  * \param[in] pme                The PME structure.
542  * \returns                      The output object.
543  */
544 GPU_FUNC_QUALIFIER PmeOutput
545 pme_gpu_getEnergyAndVirial(const gmx_pme_t &GPU_FUNC_ARGUMENT(pme)) GPU_FUNC_TERM_WITH_RETURN(PmeOutput {});
546
547 /*! \libinternal \brief
548  * Returns the GPU outputs (forces, energy and virial)
549  *
550  * \param[in] pme                The PME structure.
551  * \param[in] flags              The combination of flags that affected this PME computation.
552  *                               The flags are the GMX_PME_ flags from pme.h.
553  * \returns                      The output object.
554  */
555 GPU_FUNC_QUALIFIER PmeOutput
556     pme_gpu_getOutput(const gmx_pme_t &GPU_FUNC_ARGUMENT(pme),
557                       int              GPU_FUNC_ARGUMENT(flags)) GPU_FUNC_TERM_WITH_RETURN(PmeOutput {});
558
559 /*! \libinternal \brief
560  * Updates the unit cell parameters. Does not check if update is necessary - that is done in pme_gpu_prepare_computation().
561  *
562  * \param[in] pmeGpu         The PME GPU structure.
563  * \param[in] box            The unit cell box.
564  */
565 GPU_FUNC_QUALIFIER void pme_gpu_update_input_box(PmeGpu *GPU_FUNC_ARGUMENT(pmeGpu),
566                                                  const matrix GPU_FUNC_ARGUMENT(box)) GPU_FUNC_TERM;
567
568 /*! \libinternal \brief
569  * Finishes the PME GPU computation, waiting for the output forces and/or energy/virial to be copied to the host.
570  * If forces were computed, they will have arrived at the external host buffer provided to gather.
571  * If virial/energy were computed, they will have arrived into the internal staging buffer
572  * (even though that should have already happened before even launching the gather).
573  * Finally, cudaEvent_t based GPU timers get updated if enabled. They also need stream synchronization for correctness.
574  * Additionally, device-side buffers are cleared asynchronously for the next computation.
575  *
576  * \param[in] pmeGpu         The PME GPU structure.
577  */
578 void pme_gpu_finish_computation(const PmeGpu *pmeGpu);
579
580 //! A binary enum for spline data layout transformation
581 enum class PmeLayoutTransform
582 {
583     GpuToHost,
584     HostToGpu
585 };
586
587 /*! \libinternal \brief
588  * Rearranges the atom spline data between the GPU and host layouts.
589  * Only used for test purposes so far, likely to be horribly slow.
590  *
591  * \param[in]  pmeGpu     The PME GPU structure.
592  * \param[out] atc        The PME CPU atom data structure (with a single-threaded layout).
593  * \param[in]  type       The spline data type (values or derivatives).
594  * \param[in]  dimIndex   Dimension index.
595  * \param[in]  transform  Layout transform type
596  */
597 GPU_FUNC_QUALIFIER void pme_gpu_transform_spline_atom_data(const PmeGpu *GPU_FUNC_ARGUMENT(pmeGpu),
598                                                            const PmeAtomComm *GPU_FUNC_ARGUMENT(atc),
599                                                            PmeSplineDataType GPU_FUNC_ARGUMENT(type),
600                                                            int GPU_FUNC_ARGUMENT(dimIndex),
601                                                            PmeLayoutTransform GPU_FUNC_ARGUMENT(transform)) GPU_FUNC_TERM;
602
603 /*! \libinternal \brief
604  * Gets a unique index to an element in a spline parameter buffer (theta/dtheta),
605  * which is laid out for GPU spread/gather kernels. The index is wrt the execution block,
606  * in range(0, atomsPerBlock * order * DIM).
607  * This is a wrapper, only used in unit tests.
608  * \param[in] order            PME order
609  * \param[in] splineIndex      Spline contribution index (from 0 to \p order - 1)
610  * \param[in] dimIndex         Dimension index (from 0 to 2)
611  * \param[in] atomIndex        Atom index wrt the block.
612  * \param[in] atomsPerWarp     Number of atoms processed by a warp.
613  *
614  * \returns Index into theta or dtheta array using GPU layout.
615  */
616 int getSplineParamFullIndex(int order,
617                             int splineIndex,
618                             int dimIndex,
619                             int atomIndex,
620                             int atomsPerWarp);
621
622 /*! \libinternal \brief
623  * Get the normal/padded grid dimensions of the real-space PME grid on GPU. Only used in tests.
624  *
625  * \param[in] pmeGpu             The PME GPU structure.
626  * \param[out] gridSize          Pointer to the grid dimensions to fill in.
627  * \param[out] paddedGridSize    Pointer to the padded grid dimensions to fill in.
628  */
629 GPU_FUNC_QUALIFIER void pme_gpu_get_real_grid_sizes(const PmeGpu *GPU_FUNC_ARGUMENT(pmeGpu),
630                                                     gmx::IVec *GPU_FUNC_ARGUMENT(gridSize),
631                                                     gmx::IVec *GPU_FUNC_ARGUMENT(paddedGridSize)) GPU_FUNC_TERM;
632
633 /*! \libinternal \brief
634  * (Re-)initializes the PME GPU data at the beginning of the run or on DLB.
635  *
636  * \param[in,out] pme             The PME structure.
637  * \param[in]     gpuInfo         The GPU information structure.
638  * \param[in]     pmeGpuProgram   The PME GPU program data
639  * \throws gmx::NotImplementedError if this generally valid PME structure is not valid for GPU runs.
640  */
641 GPU_FUNC_QUALIFIER void pme_gpu_reinit(gmx_pme_t *GPU_FUNC_ARGUMENT(pme),
642                                        const gmx_device_info_t *GPU_FUNC_ARGUMENT(gpuInfo),
643                                        PmeGpuProgramHandle GPU_FUNC_ARGUMENT(pmeGpuProgram)) GPU_FUNC_TERM;
644
645 /*! \libinternal \brief
646  * Destroys the PME GPU data at the end of the run.
647  *
648  * \param[in] pmeGpu     The PME GPU structure.
649  */
650 GPU_FUNC_QUALIFIER void pme_gpu_destroy(PmeGpu *GPU_FUNC_ARGUMENT(pmeGpu)) GPU_FUNC_TERM;
651
652 /*! \libinternal \brief
653  * Reallocates the local atoms data (charges, coordinates, etc.). Copies the charges to the GPU.
654  *
655  * \param[in] pmeGpu    The PME GPU structure.
656  * \param[in] nAtoms    The number of particles.
657  * \param[in] charges   The pointer to the host-side array of particle charges.
658  *
659  * This is a function that should only be called in the beginning of the run and on domain decomposition.
660  * Should be called before the pme_gpu_set_io_ranges.
661  */
662 GPU_FUNC_QUALIFIER void pme_gpu_reinit_atoms(PmeGpu *GPU_FUNC_ARGUMENT(pmeGpu),
663                                              int         GPU_FUNC_ARGUMENT(nAtoms),
664                                              const real       *GPU_FUNC_ARGUMENT(charges)) GPU_FUNC_TERM;
665
666 /*! \brief \libinternal
667  * The PME GPU reinitialization function that is called both at the end of any PME computation and on any load balancing.
668  *
669  * This clears the device-side working buffers in preparation for new computation.
670  *
671  * \param[in] pmeGpu            The PME GPU structure.
672  */
673 void pme_gpu_reinit_computation(const PmeGpu *pmeGpu);
674
675 /*! \brief
676  * Blocks until PME GPU tasks are completed, and gets the output forces and virial/energy
677  * (if they were to be computed).
678  *
679  * \param[in]  pme            The PME data structure.
680  * \param[in]  flags          The combination of flags to affect this PME computation.
681  *                            The flags are the GMX_PME_ flags from pme.h.
682  * \param[out] wcycle         The wallclock counter.
683  * \return     The output forces, energy and virial
684  */
685 GPU_FUNC_QUALIFIER PmeOutput
686     pme_gpu_wait_finish_task(gmx_pme_t            *GPU_FUNC_ARGUMENT(pme),
687                              int                   GPU_FUNC_ARGUMENT(flags),
688                              gmx_wallcycle        *GPU_FUNC_ARGUMENT(wcycle)) GPU_FUNC_TERM_WITH_RETURN(PmeOutput {}
689                                                                                                         );
690
691 #endif