Pipeline GPU PME Spline/Spread with PP Comms
[alexxy/gromacs.git] / src / gromacs / ewald / pme.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
7  * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 /*! \libinternal \file
39  *
40  * \brief This file contains function declarations necessary for
41  * computing energies and forces for the PME long-ranged part (Coulomb
42  * and LJ).
43  *
44  * \author Berk Hess <hess@kth.se>
45  * \inlibraryapi
46  * \ingroup module_ewald
47  */
48
49 #ifndef GMX_EWALD_PME_H
50 #define GMX_EWALD_PME_H
51
52 #include <string>
53 #include <vector>
54
55 #include "gromacs/gpu_utils/devicebuffer_datatype.h"
56 #include "gromacs/gpu_utils/gpu_macros.h"
57 #include "gromacs/math/vectypes.h"
58 #include "gromacs/utility/real.h"
59
60 struct gmx_hw_info_t;
61 struct t_commrec;
62 struct t_inputrec;
63 struct t_nrnb;
64 struct PmeGpu;
65 struct gmx_wallclock_gpu_pme_t;
66 struct gmx_enerdata_t;
67 struct gmx_mtop_t;
68 struct gmx_pme_t;
69 struct gmx_wallcycle;
70 struct NumPmeDomains;
71
72 class DeviceContext;
73 class DeviceStream;
74 enum class GpuTaskCompletion;
75 class PmeGpuProgram;
76 class GpuEventSynchronizer;
77
78 /*! \brief Hack to selectively enable some parts of PME during unit testing.
79  *
80  * Set to \c false by default. If any of the tests sets it to \c true, it will
81  * make the compatibility check consider PME to be supported in SYCL builds.
82  *
83  * Currently we don't have proper PME implementation with SYCL, but we still want
84  * to run tests for some of the kernels.
85  *
86  * \todo Remove after #3927 is done and PME is fully enabled in SYCL builds.
87  */
88 //NOLINTNEXTLINE(cppcoreguidelines-avoid-non-const-global-variables)
89 extern bool g_allowPmeWithSyclForTesting;
90
91 namespace gmx
92 {
93 template<typename>
94 class ArrayRef;
95 class ForceWithVirial;
96 class MDLogger;
97 enum class PinningPolicy : int;
98 class StepWorkload;
99
100 /*! \libinternal \brief Class for managing usage of separate PME-only ranks
101  *
102  * Used for checking if some parts of the code could not use PME-only ranks
103  *
104  */
105 class SeparatePmeRanksPermitted
106 {
107 public:
108     //! Disables PME ranks permitted flag with a reason
109     void disablePmeRanks(const std::string& reason);
110     //! Return status of PME ranks usage
111     bool permitSeparatePmeRanks() const;
112     //! Returns all reasons, for not using PME ranks
113     std::string reasonsWhyDisabled() const;
114
115 private:
116     //! Flag that informs whether simualtion could use dedicated PME ranks
117     bool permitSeparatePmeRanks_ = true;
118     //! Storage for all reasons, why PME ranks could not be used
119     std::vector<std::string> reasons_;
120 };
121
122 class PmeCoordinateReceiverGpu;
123 } // namespace gmx
124
125 enum
126 {
127     GMX_SUM_GRID_FORWARD,
128     GMX_SUM_GRID_BACKWARD
129 };
130
131 /*! \brief Possible PME codepaths on a rank.
132  * \todo: make this enum class with gmx_pme_t C++ refactoring
133  */
134 enum class PmeRunMode
135 {
136     None,  //!< No PME task is done
137     CPU,   //!< Whole PME computation is done on CPU
138     GPU,   //!< Whole PME computation is done on GPU
139     Mixed, //!< Mixed mode: only spread and gather run on GPU; FFT and solving are done on CPU.
140 };
141
142 /*! \brief Return the smallest allowed PME grid size for \p pmeOrder */
143 int minimalPmeGridSize(int pmeOrder);
144
145 //! Return whether the grid of \c pme is identical to \c grid_size.
146 bool gmx_pme_grid_matches(const gmx_pme_t& pme, const ivec grid_size);
147
148 /*! \brief Check restrictions on pme_order and the PME grid nkx,nky,nkz.
149  *
150  * With errorsAreFatal=true, an exception or fatal error is generated
151  * on violation of restrictions.
152  * With errorsAreFatal=false, false is returned on violation of restrictions.
153  * When all restrictions are obeyed, true is returned.
154  * Argument useThreads tells if any MPI rank doing PME uses more than 1 threads.
155  * If at calling useThreads is unknown, pass true for conservative checking.
156  *
157  * The PME GPU restrictions are checked separately during pme_gpu_init().
158  */
159 bool gmx_pme_check_restrictions(int  pme_order,
160                                 int  nkx,
161                                 int  nky,
162                                 int  nkz,
163                                 int  numPmeDomainsAlongX,
164                                 bool useThreads,
165                                 bool errorsAreFatal);
166
167 /*! \brief Construct PME data
168  *
169  * \throws   gmx::InconsistentInputError if input grid sizes/PME order are inconsistent.
170  * \returns  Pointer to newly allocated and initialized PME data.
171  *
172  * \todo We should evolve something like a \c GpuManager that holds \c
173  * DeviceInformation* and \c PmeGpuProgram* and perhaps other
174  * related things whose lifetime can/should exceed that of a task (or
175  * perhaps task manager). See Issue #2522.
176  */
177 gmx_pme_t* gmx_pme_init(const t_commrec*     cr,
178                         const NumPmeDomains& numPmeDomains,
179                         const t_inputrec*    ir,
180                         gmx_bool             bFreeEnergy_q,
181                         gmx_bool             bFreeEnergy_lj,
182                         gmx_bool             bReproducible,
183                         real                 ewaldcoeff_q,
184                         real                 ewaldcoeff_lj,
185                         int                  nthread,
186                         PmeRunMode           runMode,
187                         PmeGpu*              pmeGpu,
188                         const DeviceContext* deviceContext,
189                         const DeviceStream*  deviceStream,
190                         const PmeGpuProgram* pmeGpuProgram,
191                         const gmx::MDLogger& mdlog);
192
193 /*! \brief As gmx_pme_init, but takes most settings, except the grid/Ewald coefficients, from
194  * pme_src. This is only called when the PME cut-off/grid size changes.
195  */
196 void gmx_pme_reinit(gmx_pme_t**       pmedata,
197                     const t_commrec*  cr,
198                     gmx_pme_t*        pme_src,
199                     const t_inputrec* ir,
200                     const ivec        grid_size,
201                     real              ewaldcoeff_q,
202                     real              ewaldcoeff_lj);
203
204 /*! \brief Destroys the PME data structure.*/
205 void gmx_pme_destroy(gmx_pme_t* pme);
206
207 /*! \brief Do a PME calculation on a CPU for the long range electrostatics and/or LJ.
208  *
209  * Computes the PME forces and the energy and viral, when requested,
210  * for all atoms in \p coordinates. Forces, when requested, are added
211  * to the buffer \p forces, which is allowed to contain more elements
212  * than the number of elements in \p coordinates.
213  * The meaning of \p flags is defined above, and determines which
214  * parts of the calculation are performed.
215  *
216  * \return 0 indicates all well, non zero is an error code.
217  */
218 int gmx_pme_do(struct gmx_pme_t*              pme,
219                gmx::ArrayRef<const gmx::RVec> coordinates,
220                gmx::ArrayRef<gmx::RVec>       forces,
221                gmx::ArrayRef<const real>      chargeA,
222                gmx::ArrayRef<const real>      chargeB,
223                gmx::ArrayRef<const real>      c6A,
224                gmx::ArrayRef<const real>      c6B,
225                gmx::ArrayRef<const real>      sigmaA,
226                gmx::ArrayRef<const real>      sigmaB,
227                const matrix                   box,
228                const t_commrec*               cr,
229                int                            maxshift_x,
230                int                            maxshift_y,
231                t_nrnb*                        nrnb,
232                gmx_wallcycle*                 wcycle,
233                matrix                         vir_q,
234                matrix                         vir_lj,
235                real*                          energy_q,
236                real*                          energy_lj,
237                real                           lambda_q,
238                real                           lambda_lj,
239                real*                          dvdlambda_q,
240                real*                          dvdlambda_lj,
241                const gmx::StepWorkload&       stepWork);
242
243 /*! \brief Calculate the PME grid energy V for n charges.
244  *
245  * The potential (found in \p pme) must have been found already with a
246  * call to gmx_pme_do(). Note that the charges are not spread on the grid in the
247  * pme struct. Currently does not work in parallel or with free
248  * energy.
249  */
250 real gmx_pme_calc_energy(gmx_pme_t* pme, gmx::ArrayRef<const gmx::RVec> x, gmx::ArrayRef<const real> q);
251
252 /*! \brief
253  * This function updates the local atom data on GPU after DD (charges, coordinates, etc.).
254  * TODO: it should update the PME CPU atom data as well.
255  * (currently PME CPU call gmx_pme_do() gets passed the input pointers for each computation).
256  *
257  * \param[in,out] pme        The PME structure.
258  * \param[in]     numAtoms   The number of particles.
259  * \param[in]     chargesA   The pointer to the array of particle charges in the normal state or FEP
260  * state A. Can be nullptr if PME is not performed on the GPU.
261  * \param[in]     chargesB   The pointer to the array of particle charges in state B. Only used if
262  * charges are perturbed and can otherwise be nullptr.
263  */
264 void gmx_pme_reinit_atoms(gmx_pme_t*                pme,
265                           int                       numAtoms,
266                           gmx::ArrayRef<const real> chargesA,
267                           gmx::ArrayRef<const real> chargesB);
268
269 /* A block of PME GPU functions */
270
271 /*! \brief Checks whether the GROMACS build allows to run PME on GPU.
272  * TODO: this partly duplicates an internal PME assert function
273  * pme_gpu_check_restrictions(), except that works with a
274  * formed gmx_pme_t structure. Should that one go away/work with inputrec?
275  *
276  * \param[out] error   If non-null, the error message when PME is not supported on GPU.
277  *
278  * \returns true if PME can run on GPU on this build, false otherwise.
279  */
280 bool pme_gpu_supports_build(std::string* error);
281
282 /*! \brief Checks whether the detected (GPU) hardware allows to run PME on GPU.
283  *
284  * \param[in]  hwinfo  Information about the detected hardware
285  * \param[out] error   If non-null, the error message when PME is not supported on GPU.
286  *
287  * \returns true if PME can run on GPU on this build, false otherwise.
288  */
289 bool pme_gpu_supports_hardware(const gmx_hw_info_t& hwinfo, std::string* error);
290
291 /*! \brief Checks whether the input system allows to run PME on GPU.
292  * TODO: this partly duplicates an internal PME assert function
293  * pme_gpu_check_restrictions(), except that works with a
294  * formed gmx_pme_t structure. Should that one go away/work with inputrec?
295  *
296  * \param[in]  ir     Input system.
297  * \param[out] error  If non-null, the error message if the input is not supported on GPU.
298  *
299  * \returns true if PME can run on GPU with this input, false otherwise.
300  */
301 bool pme_gpu_supports_input(const t_inputrec& ir, std::string* error);
302
303 /*! \brief
304  * Returns the active PME codepath (CPU, GPU, mixed).
305  * \todo This is a rather static data that should be managed by the higher level task scheduler.
306  *
307  * \param[in]  pme            The PME data structure.
308  * \returns active PME codepath.
309  */
310 PmeRunMode pme_run_mode(const gmx_pme_t* pme);
311
312 /*! \libinternal \brief
313  * Return the pinning policy appropriate for this build configuration
314  * for relevant buffers used for PME task on this rank (e.g. running
315  * on a GPU). */
316 gmx::PinningPolicy pme_get_pinning_policy();
317
318 /*! \brief
319  * Tells if PME is enabled to run on GPU (not necessarily active at the moment).
320  * \todo This is a rather static data that should be managed by the hardware assignment manager.
321  * For now, it is synonymous with the active PME codepath (in the absence of dynamic switching).
322  *
323  * \param[in]  pme            The PME data structure.
324  * \returns true if PME can run on GPU, false otherwise.
325  */
326 inline bool pme_gpu_task_enabled(const gmx_pme_t* pme)
327 {
328     return (pme != nullptr) && (pme_run_mode(pme) != PmeRunMode::CPU);
329 }
330
331 /*! \brief Returns the block size requirement
332  *
333  * The GPU version of PME requires that the coordinates array have a
334  * size divisible by the returned number.
335  *
336  * \param[in]  pme  The PME data structure.
337  */
338 GPU_FUNC_QUALIFIER int pme_gpu_get_block_size(const gmx_pme_t* GPU_FUNC_ARGUMENT(pme))
339         GPU_FUNC_TERM_WITH_RETURN(0);
340
341 // The following functions are all the PME GPU entry points,
342 // currently inlining to nothing on non-CUDA builds.
343
344 /*! \brief
345  * Resets the PME GPU timings. To be called at the reset step.
346  *
347  * \param[in] pme            The PME structure.
348  */
349 GPU_FUNC_QUALIFIER void pme_gpu_reset_timings(const gmx_pme_t* GPU_FUNC_ARGUMENT(pme)) GPU_FUNC_TERM;
350
351 /*! \brief
352  * Copies the PME GPU timings to the gmx_wallclock_gpu_pme_t structure (for log output). To be called at the run end.
353  *
354  * \param[in] pme               The PME structure.
355  * \param[in] timings           The gmx_wallclock_gpu_pme_t structure.
356  */
357 GPU_FUNC_QUALIFIER void pme_gpu_get_timings(const gmx_pme_t* GPU_FUNC_ARGUMENT(pme),
358                                             gmx_wallclock_gpu_pme_t* GPU_FUNC_ARGUMENT(timings)) GPU_FUNC_TERM;
359
360 /* The main PME GPU functions */
361
362 /*! \brief
363  * Prepares PME on GPU computation (updating the box if needed)
364  * \param[in] pme               The PME data structure.
365  * \param[in] box               The unit cell box.
366  * \param[in] wcycle            The wallclock counter.
367  * \param[in] stepWork          The required work for this simulation step
368  */
369 GPU_FUNC_QUALIFIER void pme_gpu_prepare_computation(gmx_pme_t*     GPU_FUNC_ARGUMENT(pme),
370                                                     const matrix   GPU_FUNC_ARGUMENT(box),
371                                                     gmx_wallcycle* GPU_FUNC_ARGUMENT(wcycle),
372                                                     const gmx::StepWorkload& GPU_FUNC_ARGUMENT(stepWork)) GPU_FUNC_TERM;
373
374 /*! \brief
375  * Launches first stage of PME on GPU - spreading kernel.
376  *
377  * \param[in] pme                            The PME data structure.
378  * \param[in] xReadyOnDevice                 Event synchronizer indicating that the coordinates
379  *                                           are ready in the device memory; nullptr allowed only
380  *                                           on separate PME ranks.
381  * \param[in] wcycle                         The wallclock counter.
382  * \param[in] lambdaQ                        The Coulomb lambda of the current state of the
383  *                                           system. Only used if FEP of Coulomb is active.
384  * \param[in] useGpuDirectComm               Whether direct GPU PME-PP communication is active
385  * \param[in]  pmeCoordinateReceiverGpu      Coordinate receiver object, which must be valid when
386  *                                           direct GPU PME-PP communication is active
387  */
388 GPU_FUNC_QUALIFIER void pme_gpu_launch_spread(
389         gmx_pme_t*                     GPU_FUNC_ARGUMENT(pme),
390         GpuEventSynchronizer*          GPU_FUNC_ARGUMENT(xReadyOnDevice),
391         gmx_wallcycle*                 GPU_FUNC_ARGUMENT(wcycle),
392         real                           GPU_FUNC_ARGUMENT(lambdaQ),
393         const bool                     GPU_FUNC_ARGUMENT(useGpuDirectComm),
394         gmx::PmeCoordinateReceiverGpu* GPU_FUNC_ARGUMENT(pmeCoordinateReceiverGpu)) GPU_FUNC_TERM;
395
396 /*! \brief
397  * Launches middle stages of PME (FFT R2C, solving, FFT C2R) either on GPU or on CPU, depending on the run mode.
398  *
399  * \param[in] pme               The PME data structure.
400  * \param[in] wcycle            The wallclock counter.
401  * \param[in] stepWork          The required work for this simulation step
402  */
403 GPU_FUNC_QUALIFIER void
404 pme_gpu_launch_complex_transforms(gmx_pme_t*               GPU_FUNC_ARGUMENT(pme),
405                                   gmx_wallcycle*           GPU_FUNC_ARGUMENT(wcycle),
406                                   const gmx::StepWorkload& GPU_FUNC_ARGUMENT(stepWork)) GPU_FUNC_TERM;
407
408 /*! \brief
409  * Launches last stage of PME on GPU - force gathering and D2H force transfer.
410  *
411  * \param[in] pme               The PME data structure.
412  * \param[in] wcycle            The wallclock counter.
413  * \param[in] lambdaQ           The Coulomb lambda to use when calculating the results.
414  */
415 GPU_FUNC_QUALIFIER void pme_gpu_launch_gather(const gmx_pme_t* GPU_FUNC_ARGUMENT(pme),
416                                               gmx_wallcycle*   GPU_FUNC_ARGUMENT(wcycle),
417                                               real GPU_FUNC_ARGUMENT(lambdaQ)) GPU_FUNC_TERM;
418
419 /*! \brief
420  * Attempts to complete PME GPU tasks.
421  *
422  * The \p completionKind argument controls whether the function blocks until all
423  * PME GPU tasks enqueued completed (as pme_gpu_wait_finish_task() does) or only
424  * checks and returns immediately if they did not.
425  * When blocking or the tasks have completed it also gets the output forces
426  * by assigning the ArrayRef to the \p forces pointer passed in.
427  * Virial/energy are also outputs if they were to be computed.
428  *
429  * \param[in]  pme             The PME data structure.
430  * \param[in]  stepWork        The required work for this simulation step
431  * \param[in]  wcycle          The wallclock counter.
432  * \param[out] forceWithVirial The output force and virial
433  * \param[out] enerd           The output energies
434  * \param[in]  lambdaQ         The Coulomb lambda to use when calculating the results.
435  * \param[in]  completionKind  Indicates whether PME task completion should only be checked rather
436  *                             than waited for
437  * \returns                    True if the PME GPU tasks have completed
438  */
439 GPU_FUNC_QUALIFIER bool pme_gpu_try_finish_task(gmx_pme_t*               GPU_FUNC_ARGUMENT(pme),
440                                                 const gmx::StepWorkload& GPU_FUNC_ARGUMENT(stepWork),
441                                                 gmx_wallcycle*           GPU_FUNC_ARGUMENT(wcycle),
442                                                 gmx::ForceWithVirial* GPU_FUNC_ARGUMENT(forceWithVirial),
443                                                 gmx_enerdata_t*       GPU_FUNC_ARGUMENT(enerd),
444                                                 real                  GPU_FUNC_ARGUMENT(lambdaQ),
445                                                 GpuTaskCompletion GPU_FUNC_ARGUMENT(completionKind))
446         GPU_FUNC_TERM_WITH_RETURN(false);
447
448 /*! \brief
449  * Blocks until PME GPU tasks are completed, and gets the output forces and virial/energy
450  * (if they were to be computed).
451  *
452  * \param[in]  pme             The PME data structure.
453  * \param[in]  stepWork        The required work for this simulation step
454  * \param[in]  wcycle          The wallclock counter.
455  * \param[out] forceWithVirial The output force and virial
456  * \param[out] enerd           The output energies
457  * \param[in]  lambdaQ         The Coulomb lambda to use when calculating the results.
458  */
459 GPU_FUNC_QUALIFIER void pme_gpu_wait_and_reduce(gmx_pme_t*               GPU_FUNC_ARGUMENT(pme),
460                                                 const gmx::StepWorkload& GPU_FUNC_ARGUMENT(stepWork),
461                                                 gmx_wallcycle*           GPU_FUNC_ARGUMENT(wcycle),
462                                                 gmx::ForceWithVirial* GPU_FUNC_ARGUMENT(forceWithVirial),
463                                                 gmx_enerdata_t*       GPU_FUNC_ARGUMENT(enerd),
464                                                 real GPU_FUNC_ARGUMENT(lambdaQ)) GPU_FUNC_TERM;
465
466 /*! \brief
467  * The PME GPU reinitialization function that is called both at the end of any PME computation and on any load balancing.
468  *
469  * Clears the internal grid and energy/virial buffers; it is not safe to start
470  * the PME computation without calling this.
471  * Note that unlike in the nbnxn module, the force buffer does not need clearing.
472  *
473  * \todo Rename this function to *clear* -- it clearly only does output resetting
474  * and we should be clear about what the function does..
475  *
476  * \param[in] pme            The PME data structure.
477  * \param[in] wcycle         The wallclock counter.
478  */
479 GPU_FUNC_QUALIFIER void pme_gpu_reinit_computation(const gmx_pme_t* GPU_FUNC_ARGUMENT(pme),
480                                                    gmx_wallcycle* GPU_FUNC_ARGUMENT(wcycle)) GPU_FUNC_TERM;
481
482 /*! \brief Set pointer to device copy of coordinate data.
483  * \param[in] pme            The PME data structure.
484  * \param[in] d_x            The pointer to the positions buffer to be set
485  */
486 GPU_FUNC_QUALIFIER void pme_gpu_set_device_x(const gmx_pme_t*        GPU_FUNC_ARGUMENT(pme),
487                                              DeviceBuffer<gmx::RVec> GPU_FUNC_ARGUMENT(d_x)) GPU_FUNC_TERM;
488
489 /*! \brief Get pointer to device copy of force data.
490  * \param[in] pme            The PME data structure.
491  * \returns                  Pointer to force data
492  */
493 GPU_FUNC_QUALIFIER DeviceBuffer<gmx::RVec> pme_gpu_get_device_f(const gmx_pme_t* GPU_FUNC_ARGUMENT(pme))
494         GPU_FUNC_TERM_WITH_RETURN(DeviceBuffer<gmx::RVec>{});
495
496 /*! \brief Get pointer to the device synchronizer object that allows syncing on PME force calculation completion
497  * \param[in] pme            The PME data structure.
498  * \returns                  Pointer to synchronizer
499  */
500 GPU_FUNC_QUALIFIER GpuEventSynchronizer* pme_gpu_get_f_ready_synchronizer(const gmx_pme_t* GPU_FUNC_ARGUMENT(pme))
501         GPU_FUNC_TERM_WITH_RETURN(nullptr);
502
503 #endif