Merge release-2019 into master
[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, 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/gpu_macros.h"
54 #include "gromacs/math/vectypes.h"
55 #include "gromacs/timing/walltime_accounting.h"
56 #include "gromacs/utility/arrayref.h"
57 #include "gromacs/utility/basedefinitions.h"
58 #include "gromacs/utility/real.h"
59
60 struct gmx_hw_info_t;
61 struct interaction_const_t;
62 struct t_commrec;
63 struct t_inputrec;
64 struct t_nrnb;
65 struct PmeGpu;
66 struct gmx_wallclock_gpu_pme_t;
67 struct gmx_device_info_t;
68 struct gmx_mtop_t;
69 struct gmx_pme_t;
70 struct gmx_wallcycle;
71 struct NumPmeDomains;
72
73 enum class GpuTaskCompletion;
74 class PmeGpuProgram;
75 //! Convenience name.
76 using PmeGpuProgramHandle = const PmeGpuProgram *;
77
78 namespace gmx
79 {
80 class ForceWithVirial;
81 class MDLogger;
82 enum class PinningPolicy : int;
83 }
84
85 enum {
86     GMX_SUM_GRID_FORWARD, GMX_SUM_GRID_BACKWARD
87 };
88
89 /*! \brief Possible PME codepaths on a rank.
90  * \todo: make this enum class with gmx_pme_t C++ refactoring
91  */
92 enum class PmeRunMode
93 {
94     None,    //!< No PME task is done
95     CPU,     //!< Whole PME computation is done on CPU
96     GPU,     //!< Whole PME computation is done on GPU
97     Mixed,   //!< Mixed mode: only spread and gather run on GPU; FFT and solving are done on CPU.
98 };
99
100 //! PME gathering output forces treatment
101 enum class PmeForceOutputHandling
102 {
103     Set,             /**< Gather simply writes into provided force buffer */
104     ReduceWithInput, /**< Gather adds its output to the buffer.
105                         On GPU, that means additional H2D copy before the kernel launch. */
106 };
107
108 /*! \brief Return the smallest allowed PME grid size for \p pmeOrder */
109 int minimalPmeGridSize(int pmeOrder);
110
111 /*! \brief Check restrictions on pme_order and the PME grid nkx,nky,nkz.
112  *
113  * With errorsAreFatal=true, an exception or fatal error is generated
114  * on violation of restrictions.
115  * With errorsAreFatal=false, false is returned on violation of restrictions.
116  * When all restrictions are obeyed, true is returned.
117  * Argument useThreads tells if any MPI rank doing PME uses more than 1 threads.
118  * If at calling useThreads is unknown, pass true for conservative checking.
119  *
120  * The PME GPU restrictions are checked separately during pme_gpu_init().
121  */
122 bool gmx_pme_check_restrictions(int pme_order,
123                                 int nkx, int nky, int nkz,
124                                 int numPmeDomainsAlongX,
125                                 bool useThreads,
126                                 bool errorsAreFatal);
127
128 /*! \brief Construct PME data
129  *
130  * \throws   gmx::InconsistentInputError if input grid sizes/PME order are inconsistent.
131  * \returns  Pointer to newly allocated and initialized PME data.
132  *
133  * \todo We should evolve something like a \c GpuManager that holds \c
134  * gmx_device_info_t * and \c PmeGpuProgramHandle and perhaps other
135  * related things whose lifetime can/should exceed that of a task (or
136  * perhaps task manager). See Redmine #2522.
137  */
138 gmx_pme_t *gmx_pme_init(const t_commrec *cr,
139                         const NumPmeDomains &numPmeDomains,
140                         const t_inputrec *ir, int homenr,
141                         gmx_bool bFreeEnergy_q, gmx_bool bFreeEnergy_lj,
142                         gmx_bool bReproducible,
143                         real ewaldcoeff_q, real ewaldcoeff_lj,
144                         int nthread,
145                         PmeRunMode runMode,
146                         PmeGpu *pmeGpu,
147                         const gmx_device_info_t *gpuInfo,
148                         PmeGpuProgramHandle pmeGpuProgram,
149                         const gmx::MDLogger &mdlog);
150
151 /*! \brief Destroys the PME data structure.*/
152 void gmx_pme_destroy(gmx_pme_t *pme);
153
154 //@{
155 /*! \brief Flag values that control what gmx_pme_do() will calculate
156  *
157  * These can be combined with bitwise-OR if more than one thing is required.
158  */
159 #define GMX_PME_SPREAD        (1<<0)
160 #define GMX_PME_SOLVE         (1<<1)
161 #define GMX_PME_CALC_F        (1<<2)
162 #define GMX_PME_CALC_ENER_VIR (1<<3)
163 /* This forces the grid to be backtransformed even without GMX_PME_CALC_F */
164 #define GMX_PME_CALC_POT      (1<<4)
165
166 #define GMX_PME_DO_ALL_F  (GMX_PME_SPREAD | GMX_PME_SOLVE | GMX_PME_CALC_F)
167 //@}
168
169 /*! \brief Do a PME calculation on a CPU for the long range electrostatics and/or LJ.
170  *
171  * The meaning of \p flags is defined above, and determines which
172  * parts of the calculation are performed.
173  *
174  * \return 0 indicates all well, non zero is an error code.
175  */
176 int gmx_pme_do(struct gmx_pme_t *pme,
177                int start,       int homenr,
178                rvec x[],        rvec f[],
179                real chargeA[],  real chargeB[],
180                real c6A[],      real c6B[],
181                real sigmaA[],   real sigmaB[],
182                matrix box,      const t_commrec *cr,
183                int  maxshift_x, int maxshift_y,
184                t_nrnb *nrnb,    gmx_wallcycle *wcycle,
185                matrix vir_q,    matrix vir_lj,
186                real *energy_q,  real *energy_lj,
187                real lambda_q,   real lambda_lj,
188                real *dvdlambda_q, real *dvdlambda_lj,
189                int flags);
190
191 /*! \brief Called on the nodes that do PME exclusively */
192 int gmx_pmeonly(struct gmx_pme_t *pme,
193                 const t_commrec *cr,     t_nrnb *mynrnb,
194                 gmx_wallcycle  *wcycle,
195                 gmx_walltime_accounting_t walltime_accounting,
196                 t_inputrec *ir, PmeRunMode runMode);
197
198 /*! \brief Calculate the PME grid energy V for n charges.
199  *
200  * The potential (found in \p pme) must have been found already with a
201  * call to gmx_pme_do() with at least GMX_PME_SPREAD and GMX_PME_SOLVE
202  * specified. Note that the charges are not spread on the grid in the
203  * pme struct. Currently does not work in parallel or with free
204  * energy.
205  */
206 void gmx_pme_calc_energy(struct gmx_pme_t *pme, int n, rvec *x, real *q, real *V);
207
208 /*! \brief Send the charges and maxshift to out PME-only node. */
209 void gmx_pme_send_parameters(const t_commrec *cr,
210                              const interaction_const_t *ic,
211                              gmx_bool bFreeEnergy_q, gmx_bool bFreeEnergy_lj,
212                              real *chargeA, real *chargeB,
213                              real *sqrt_c6A, real *sqrt_c6B,
214                              real *sigmaA, real *sigmaB,
215                              int maxshift_x, int maxshift_y);
216
217 /*! \brief Send the coordinates to our PME-only node and request a PME calculation */
218 void gmx_pme_send_coordinates(const t_commrec *cr, matrix box, rvec *x,
219                               real lambda_q, real lambda_lj,
220                               gmx_bool bEnerVir,
221                               int64_t step, gmx_wallcycle *wcycle);
222
223 /*! \brief Tell our PME-only node to finish */
224 void gmx_pme_send_finish(const t_commrec *cr);
225
226 /*! \brief Tell our PME-only node to reset all cycle and flop counters */
227 void gmx_pme_send_resetcounters(const t_commrec *cr, int64_t step);
228
229 /*! \brief PP nodes receive the long range forces from the PME nodes */
230 void gmx_pme_receive_f(const t_commrec *cr,
231                        gmx::ForceWithVirial *forceWithVirial,
232                        real *energy_q, real *energy_lj,
233                        real *dvdlambda_q, real *dvdlambda_lj,
234                        float *pme_cycles);
235
236 /*! \brief
237  * This function updates the local atom data on GPU after DD (charges, coordinates, etc.).
238  * TODO: it should update the PME CPU atom data as well.
239  * (currently PME CPU call gmx_pme_do() gets passed the input pointers for each computation).
240  *
241  * \param[in] pme            The PME structure.
242  * \param[in] nAtoms         The number of particles.
243  * \param[in] charges        The pointer to the array of particle charges.
244  */
245 void gmx_pme_reinit_atoms(const gmx_pme_t *pme, int nAtoms, const real *charges);
246
247 /* A block of PME GPU functions */
248
249 /*! \brief Checks whether the GROMACS build allows to run PME on GPU.
250  * TODO: this partly duplicates an internal PME assert function
251  * pme_gpu_check_restrictions(), except that works with a
252  * formed gmx_pme_t structure. Should that one go away/work with inputrec?
253  *
254  * \param[in]  hwinfo  Information about the detected hardware
255  * \param[out] error   If non-null, the error message when PME is not supported on GPU.
256  *
257  * \returns true if PME can run on GPU on this build, false otherwise.
258  */
259 bool pme_gpu_supports_build(const gmx_hw_info_t &hwinfo,
260                             std::string         *error);
261
262 /*! \brief Checks whether the input system allows to run PME on GPU.
263  * TODO: this partly duplicates an internal PME assert function
264  * pme_gpu_check_restrictions(), except that works with a
265  * formed gmx_pme_t structure. Should that one go away/work with inputrec?
266  *
267  * \param[in]  ir     Input system.
268  * \param[in]  mtop   Complete system topology to check if an FE simulation perturbs charges.
269  * \param[out] error  If non-null, the error message if the input is not supported on GPU.
270  *
271  * \returns true if PME can run on GPU with this input, false otherwise.
272  */
273 bool pme_gpu_supports_input(const t_inputrec &ir, const gmx_mtop_t &mtop, std::string *error);
274
275 /*! \brief
276  * Returns the active PME codepath (CPU, GPU, mixed).
277  * \todo This is a rather static data that should be managed by the higher level task scheduler.
278  *
279  * \param[in]  pme            The PME data structure.
280  * \returns active PME codepath.
281  */
282 PmeRunMode pme_run_mode(const gmx_pme_t *pme);
283
284 /*! \libinternal \brief
285  * Return the pinning policy appropriate for this build configuration
286  * for relevant buffers used for PME task on this rank (e.g. running
287  * on a GPU). */
288 gmx::PinningPolicy pme_get_pinning_policy();
289
290 /*! \brief
291  * Tells if PME is enabled to run on GPU (not necessarily active at the moment).
292  * \todo This is a rather static data that should be managed by the hardware assignment manager.
293  * For now, it is synonymous with the active PME codepath (in the absence of dynamic switching).
294  *
295  * \param[in]  pme            The PME data structure.
296  * \returns true if PME can run on GPU, false otherwise.
297  */
298 inline bool pme_gpu_task_enabled(const gmx_pme_t *pme)
299 {
300     return (pme != nullptr) && (pme_run_mode(pme) != PmeRunMode::CPU);
301 }
302
303 // The following functions are all the PME GPU entry points,
304 // currently inlining to nothing on non-CUDA builds.
305
306 /*! \brief
307  * Resets the PME GPU timings. To be called at the reset step.
308  *
309  * \param[in] pme            The PME structure.
310  */
311 GPU_FUNC_QUALIFIER void pme_gpu_reset_timings(const gmx_pme_t *GPU_FUNC_ARGUMENT(pme)) GPU_FUNC_TERM
312
313 /*! \brief
314  * Copies the PME GPU timings to the gmx_wallclock_gpu_pme_t structure (for log output). To be called at the run end.
315  *
316  * \param[in] pme               The PME structure.
317  * \param[in] timings           The gmx_wallclock_gpu_pme_t structure.
318  */
319 GPU_FUNC_QUALIFIER void pme_gpu_get_timings(const gmx_pme_t         *GPU_FUNC_ARGUMENT(pme),
320                                             gmx_wallclock_gpu_pme_t *GPU_FUNC_ARGUMENT(timings)) GPU_FUNC_TERM
321
322 /* The main PME GPU functions */
323
324 /*! \brief
325  * Prepares PME on GPU computation (updating the box if needed)
326  * \param[in] pme               The PME data structure.
327  * \param[in] needToUpdateBox   Tells if the stored unit cell parameters should be updated from \p box.
328  * \param[in] box               The unit cell box.
329  * \param[in] wcycle            The wallclock counter.
330  * \param[in] flags             The combination of flags to affect this PME computation.
331  *                              The flags are the GMX_PME_ flags from pme.h.
332  */
333 GPU_FUNC_QUALIFIER void pme_gpu_prepare_computation(gmx_pme_t      *GPU_FUNC_ARGUMENT(pme),
334                                                     bool            GPU_FUNC_ARGUMENT(needToUpdateBox),
335                                                     const matrix    GPU_FUNC_ARGUMENT(box),
336                                                     gmx_wallcycle  *GPU_FUNC_ARGUMENT(wcycle),
337                                                     int             GPU_FUNC_ARGUMENT(flags)) GPU_FUNC_TERM
338
339 /*! \brief
340  * Launches first stage of PME on GPU - H2D input transfers, spreading kernel, and D2H grid transfer if needed.
341  *
342  * \param[in] pme               The PME data structure.
343  * \param[in] x                 The array of local atoms' coordinates.
344  * \param[in] wcycle            The wallclock counter.
345  */
346 GPU_FUNC_QUALIFIER void pme_gpu_launch_spread(gmx_pme_t      *GPU_FUNC_ARGUMENT(pme),
347                                               const rvec     *GPU_FUNC_ARGUMENT(x),
348                                               gmx_wallcycle  *GPU_FUNC_ARGUMENT(wcycle)) GPU_FUNC_TERM
349
350 /*! \brief
351  * Launches middle stages of PME (FFT R2C, solving, FFT C2R) either on GPU or on CPU, depending on the run mode.
352  *
353  * \param[in] pme               The PME data structure.
354  * \param[in] wcycle            The wallclock counter.
355  */
356 GPU_FUNC_QUALIFIER void pme_gpu_launch_complex_transforms(gmx_pme_t       *GPU_FUNC_ARGUMENT(pme),
357                                                           gmx_wallcycle   *GPU_FUNC_ARGUMENT(wcycle)) GPU_FUNC_TERM
358
359 /*! \brief
360  * Launches last stage of PME on GPU - force gathering and D2H force transfer.
361  *
362  * \param[in]  pme               The PME data structure.
363  * \param[in]  wcycle            The wallclock counter.
364  * \param[in]  forceTreatment    Tells how data should be treated. The gathering kernel either stores
365  *                               the output reciprocal forces into the host array, or copies its contents to the GPU first
366  *                               and accumulates. The reduction is non-atomic.
367  */
368 GPU_FUNC_QUALIFIER void pme_gpu_launch_gather(const gmx_pme_t        *GPU_FUNC_ARGUMENT(pme),
369                                               gmx_wallcycle          *GPU_FUNC_ARGUMENT(wcycle),
370                                               PmeForceOutputHandling  GPU_FUNC_ARGUMENT(forceTreatment)) GPU_FUNC_TERM
371
372 /*! \brief
373  * Blocks until PME GPU tasks are completed, and gets the output forces and virial/energy
374  * (if they were to be computed).
375  *
376  * \param[in]  pme            The PME data structure.
377  * \param[out] wcycle         The wallclock counter.
378  * \param[out] forces         The output forces.
379  * \param[out] virial         The output virial matrix.
380  * \param[out] energy         The output energy.
381  */
382 GPU_FUNC_QUALIFIER void pme_gpu_wait_finish_task(const gmx_pme_t                *GPU_FUNC_ARGUMENT(pme),
383                                                  gmx_wallcycle                  *GPU_FUNC_ARGUMENT(wcycle),
384                                                  gmx::ArrayRef<const gmx::RVec> *GPU_FUNC_ARGUMENT(forces),
385                                                  matrix                          GPU_FUNC_ARGUMENT(virial),
386                                                  real                           *GPU_FUNC_ARGUMENT(energy)) GPU_FUNC_TERM
387 /*! \brief
388  * Attempts to complete PME GPU tasks.
389  *
390  * The \p completionKind argument controls whether the function blocks until all
391  * PME GPU tasks enqueued completed (as pme_gpu_wait_finish_task() does) or only
392  * checks and returns immediately if they did not.
393  * When blocking or the tasks have completed it also gets the output forces
394  * by assigning the ArrayRef to the \p forces pointer passed in.
395  * Virial/energy are also outputs if they were to be computed.
396  *
397  * Note: also launches the reinitalization of the PME output buffers.
398  * TODO: this should be moved out to avoid miscounting its wall-time (as wait iso launch).
399  *
400  * \param[in]  pme            The PME data structure.
401  * \param[in]  wcycle         The wallclock counter.
402  * \param[out] forces         The output forces.
403  * \param[out] virial         The output virial matrix.
404  * \param[out] energy         The output energy.
405  * \param[in]  completionKind  Indicates whether PME task completion should only be checked rather than waited for
406  * \returns                   True if the PME GPU tasks have completed
407  */
408 GPU_FUNC_QUALIFIER bool pme_gpu_try_finish_task(const gmx_pme_t                *GPU_FUNC_ARGUMENT(pme),
409                                                 gmx_wallcycle                  *GPU_FUNC_ARGUMENT(wcycle),
410                                                 gmx::ArrayRef<const gmx::RVec> *GPU_FUNC_ARGUMENT(forces),
411                                                 matrix                          GPU_FUNC_ARGUMENT(virial),
412                                                 real                           *GPU_FUNC_ARGUMENT(energy),
413                                                 GpuTaskCompletion               GPU_FUNC_ARGUMENT(completionKind)) GPU_FUNC_TERM_WITH_RETURN(false)
414
415 /*! \brief
416  * The PME GPU reinitialization function that is called both at the end of any PME computation and on any load balancing.
417  *
418  * Clears the internal grid and energy/virial buffers; it is not safe to start
419  * the PME computation without calling this.
420  * Note that unlike in the nbnxn module, the force buffer does not need clearing.
421  *
422  * \todo Rename this function to *clear* -- it clearly only does output resetting
423  * and we should be clear about what the function does..
424  *
425  * \param[in] pme            The PME data structure.
426  * \param[in] wcycle         The wallclock counter.
427  */
428 GPU_FUNC_QUALIFIER void pme_gpu_reinit_computation(const gmx_pme_t *GPU_FUNC_ARGUMENT(pme),
429                                                    gmx_wallcycle   *GPU_FUNC_ARGUMENT(wcycle)) GPU_FUNC_TERM
430
431 #endif