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