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