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