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