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