7a6491a4aac4538d34e0819d022b94a12c4fdbce
[alexxy/gromacs.git] / src / gromacs / nbnxm / nbnxm.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2017,2018,2019, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35
36 // FIXME: remove the "__" prefix in front of the group def when we move the
37 //        nonbonded code into separate dir.
38
39 /*! \libinternal \defgroup __module_nbnxm Short-range non-bonded interaction module
40  * \ingroup group_mdrun
41  *
42  * \brief Computes forces and energies for short-range pair-interactions
43  * based on the Verlet algorithm. The algorithm uses pair-lists generated
44  * at fixed intervals as well as various flavors of pair interaction kernels
45  * implemented for a wide range of CPU and GPU architectures.
46  *
47  * The module includes support for flavors of Coulomb and Lennard-Jones interaction
48  * treatment implemented for a large range of SIMD instruction sets for CPU
49  * architectures as well as in CUDA and OpenCL for GPU architectures.
50  * Additionally there is a reference CPU non-SIMD and a reference CPU
51  * for GPU pair-list setup interaction kernel.
52  *
53  * The implementation of the kernels is based on the cluster non-bonded algorithm
54  * which in the code is referred to as the NxM algorithms ("nbnxm_" prefix);
55  * for details of the algorithm see DOI:10.1016/j.cpc.2013.06.003.
56  *
57  * Algorithmically, the non-bonded computation has two different modes:
58  * A "classical" mode: generate a list every nstlist steps containing at least
59  * all atom pairs up to a distance of rlistOuter and compute pair interactions
60  * for all pairs that are within the interaction cut-off.
61  * A "dynamic pruning" mode: generate an "outer-list" up to cut-off rlistOuter
62  * every nstlist steps and prune the outer-list using a cut-off of rlistInner
63  * every nstlistPrune steps to obtain a, smaller, "inner-list". This
64  * results in fewer interaction computations and allows for a larger nstlist.
65  * On a GPU, this dynamic pruning is performed in a rolling fashion, pruning
66  * only a sub-part of the list each (second) step. This way it can often
67  * overlap with integration and constraints on the CPU.
68  * Currently a simple heuristic determines which mode will be used.
69  *
70  * TODO: add a summary list and brief descriptions of the different submodules:
71  * search, CPU kernels, GPU glue code + kernels.
72  *
73  * \author Berk Hess <hess@kth.se>
74  * \author Szilárd Páll <pall.szilard@gmail.com>
75  * \author Mark Abraham <mark.j.abraham@gmail.com>
76  * \author Anca Hamuraru <anca@streamcomputing.eu>
77  * \author Teemu Virolainen <teemu@streamcomputing.eu>
78  * \author Dimitrios Karkoulis <dimitris.karkoulis@gmail.com>
79  *
80  * TODO: add more authors!
81  */
82
83 /*! \libinternal \file
84  *
85  * \brief This file contains the public interface of the nbnxm module
86  * that implements the NxM atom cluster non-bonded algorithm to efficiently
87  * compute pair forces.
88  *
89  *
90  * \author Berk Hess <hess@kth.se>
91  * \author Szilárd Páll <pall.szilard@gmail.com>
92  *
93  * \inlibraryapi
94  * \ingroup __module_nbnxm
95  */
96
97
98 #ifndef GMX_NBNXM_NBNXM_H
99 #define GMX_NBNXM_NBNXM_H
100
101 #include <memory>
102
103 #include "gromacs/math/vectypes.h"
104 #include "gromacs/utility/arrayref.h"
105 #include "gromacs/utility/enumerationhelpers.h"
106 #include "gromacs/utility/real.h"
107
108 #include "locality.h"
109
110 // TODO: Remove this include and the two nbnxm includes above
111 #include "nbnxm_gpu.h"
112
113 struct gmx_device_info_t;
114 struct gmx_domdec_zones_t;
115 struct gmx_enerdata_t;
116 struct gmx_hw_info_t;
117 struct gmx_mtop_t;
118 struct gmx_wallcycle;
119 struct interaction_const_t;
120 struct nonbonded_verlet_t;
121 class PairSearch;
122 class PairlistSets;
123 struct t_blocka;
124 struct t_commrec;
125 struct t_lambda;
126 struct t_mdatoms;
127 struct t_nrnb;
128 struct t_forcerec;
129 struct t_inputrec;
130
131 /*! \brief Switch for whether to use GPU for buffer ops*/
132 enum class BufferOpsUseGpu
133 {
134     True,
135     False
136 };
137
138 /*! \brief Switch for whether forces should accumulate in GPU buffer ops */
139 enum class GpuBufferOpsAccumulateForce
140 {
141     True,  // Force should be accumulated and format converted
142     False, // Force should be not accumulated, just format converted
143     Null   // GPU buffer ops are not in use, so this object is not applicable
144 };
145
146
147 namespace gmx
148 {
149 class MDLogger;
150 class UpdateGroupsCog;
151 }
152
153 namespace Nbnxm
154 {
155 enum class KernelType;
156 }
157
158 namespace Nbnxm
159 {
160
161 /*! \brief Nonbonded NxN kernel types: plain C, CPU SIMD, GPU, GPU emulation */
162 enum class KernelType : int
163 {
164     NotSet = 0,
165     Cpu4x4_PlainC,
166     Cpu4xN_Simd_4xN,
167     Cpu4xN_Simd_2xNN,
168     Gpu8x8x8,
169     Cpu8x8x8_PlainC,
170     Count
171 };
172
173 /*! \brief Ewald exclusion types */
174 enum class EwaldExclusionType : int
175 {
176     NotSet = 0,
177     Table,
178     Analytical,
179     DecidedByGpuModule
180 };
181
182 /* \brief The non-bonded setup, also affects the pairlist construction kernel */
183 struct KernelSetup
184 {
185     //! The non-bonded type, also affects the pairlist construction kernel
186     KernelType         kernelType = KernelType::NotSet;
187     //! Ewald exclusion computation handling type, currently only used for CPU
188     EwaldExclusionType ewaldExclusionType = EwaldExclusionType::NotSet;
189 };
190
191 /*! \brief Return a string identifying the kernel type.
192  *
193  * \param [in] kernelType   nonbonded kernel type, takes values from the nbnxn_kernel_type enum
194  * \returns                 a string identifying the kernel corresponding to the type passed as argument
195  */
196 const char *lookup_kernel_name(Nbnxm::KernelType kernelType);
197
198 } // namespace Nbnxm
199
200 /*! \brief Flag to tell the nonbonded kernels whether to clear the force output buffers */
201 enum {
202     enbvClearFNo, enbvClearFYes
203 };
204
205 /*! \libinternal
206  *  \brief Top-level non-bonded data structure for the Verlet-type cut-off scheme. */
207 struct nonbonded_verlet_t
208 {
209     public:
210         //! Constructs an object from its components
211         nonbonded_verlet_t(std::unique_ptr<PairlistSets>      pairlistSets,
212                            std::unique_ptr<PairSearch>        pairSearch,
213                            std::unique_ptr<nbnxn_atomdata_t>  nbat,
214                            const Nbnxm::KernelSetup          &kernelSetup,
215                            gmx_nbnxn_gpu_t                   *gpu_nbv);
216
217         ~nonbonded_verlet_t();
218
219         //! Returns whether a GPU is use for the non-bonded calculations
220         bool useGpu() const
221         {
222             return kernelSetup_.kernelType == Nbnxm::KernelType::Gpu8x8x8;
223         }
224
225         //! Returns whether a GPU is emulated for the non-bonded calculations
226         bool emulateGpu() const
227         {
228             return kernelSetup_.kernelType == Nbnxm::KernelType::Cpu8x8x8_PlainC;
229         }
230
231         //! Return whether the pairlist is of simple, CPU type
232         bool pairlistIsSimple() const
233         {
234             return !useGpu() && !emulateGpu();
235         }
236
237         //! Initialize the pair list sets, TODO this should be private
238         void initPairlistSets(bool haveMultipleDomains);
239
240         //! Returns the order of the local atoms on the grid
241         gmx::ArrayRef<const int> getLocalAtomOrder() const;
242
243         //! Sets the order of the local atoms to the order grid atom ordering
244         void setLocalAtomOrder();
245
246         //! Returns the index position of the atoms on the search grid
247         gmx::ArrayRef<const int> getGridIndices() const;
248
249         //! Constructs the pairlist for the given locality
250         void constructPairlist(Nbnxm::InteractionLocality  iLocality,
251                                const t_blocka             *excl,
252                                int64_t                     step,
253                                t_nrnb                     *nrnb);
254
255         //! Updates all the atom properties in Nbnxm
256         void setAtomProperties(const t_mdatoms          &mdatoms,
257                                gmx::ArrayRef<const int>  atomInfo);
258
259         //! Updates the coordinates in Nbnxm for the given locality
260         void setCoordinates(Nbnxm::AtomLocality             locality,
261                             bool                            fillLocal,
262                             gmx::ArrayRef<const gmx::RVec>  x,
263                             BufferOpsUseGpu                 useGpu,
264                             void                           *xPmeDevicePtr,
265                             gmx_wallcycle                  *wcycle);
266
267         //! Init for GPU version of setup coordinates in Nbnxm
268         void atomdata_init_copy_x_to_nbat_x_gpu();
269
270         //! Sync the nonlocal GPU stream with dependent tasks in the local queue.
271         void insertNonlocalGpuDependency(Nbnxm::InteractionLocality interactionLocality);
272
273         //! Returns a reference to the pairlist sets
274         const PairlistSets &pairlistSets() const
275         {
276             return *pairlistSets_;
277         }
278
279         //! Returns whether step is a dynamic list pruning step, for CPU lists
280         bool isDynamicPruningStepCpu(int64_t step) const;
281
282         //! Returns whether step is a dynamic list pruning step, for GPU lists
283         bool isDynamicPruningStepGpu(int64_t step) const;
284
285         //! Dispatches the dynamic pruning kernel for the given locality, for CPU lists
286         void dispatchPruneKernelCpu(Nbnxm::InteractionLocality  iLocality,
287                                     const rvec                 *shift_vec);
288
289         //! Dispatches the dynamic pruning kernel for GPU lists
290         void dispatchPruneKernelGpu(int64_t step);
291
292         //! \brief Executes the non-bonded kernel of the GPU or launches it on the GPU
293         void dispatchNonbondedKernel(Nbnxm::InteractionLocality  iLocality,
294                                      const interaction_const_t  &ic,
295                                      int                         forceFlags,
296                                      int                         clearF,
297                                      const t_forcerec           &fr,
298                                      gmx_enerdata_t             *enerd,
299                                      t_nrnb                     *nrnb,
300                                      gmx_wallcycle              *wcycle);
301
302         //! Executes the non-bonded free-energy kernel, always runs on the CPU
303         void dispatchFreeEnergyKernel(Nbnxm::InteractionLocality  iLocality,
304                                       t_forcerec                 *fr,
305                                       rvec                        x[],
306                                       rvec                        f[],
307                                       const t_mdatoms            &mdatoms,
308                                       t_lambda                   *fepvals,
309                                       real                       *lambda,
310                                       gmx_enerdata_t             *enerd,
311                                       int                         forceFlags,
312                                       t_nrnb                     *nrnb,
313                                       gmx_wallcycle              *wcycle);
314
315         //! Add the forces stored in nbat to f, zeros the forces in nbat */
316         void atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality                 locality,
317                                       rvec                               *f,
318                                       BufferOpsUseGpu                     useGpu,
319                                       GpuBufferOpsAccumulateForce         accumulateForce,
320                                       gmx_wallcycle                      *wcycle);
321
322         /*! \brief Outer body of function to perform initialization for F buffer operations on GPU. */
323         void atomdata_init_add_nbat_f_to_f_gpu(gmx_wallcycle *wcycle);
324
325         /*! \brief H2D transfer of force buffer*/
326         void launch_copy_f_to_gpu(rvec *f, Nbnxm::AtomLocality locality);
327
328         /*! \brief D2H transfer of force buffer*/
329         void launch_copy_f_from_gpu(rvec *f, Nbnxm::AtomLocality locality);
330
331         /*! \brief Wait for GPU force reduction task and D2H transfer of its results to complete
332          *
333          * FIXME: need more details: when should be called / after which operation, etc.
334          */
335         void wait_for_gpu_force_reduction(Nbnxm::AtomLocality locality);
336
337         //! Return the kernel setup
338         const Nbnxm::KernelSetup &kernelSetup() const
339         {
340             return kernelSetup_;
341         }
342
343         //! Returns the outer radius for the pair list
344         real pairlistInnerRadius() const;
345
346         //! Returns the outer radius for the pair list
347         real pairlistOuterRadius() const;
348
349         //! Changes the pair-list outer and inner radius
350         void changePairlistRadii(real rlistOuter,
351                                  real rlistInner);
352
353         //! Set up internal flags that indicate what type of short-range work there is.
354         void setupGpuShortRangeWork(const gmx::GpuBonded             *gpuBonded,
355                                     const Nbnxm::InteractionLocality  iLocality)
356         {
357             if (useGpu() && !emulateGpu())
358             {
359                 Nbnxm::setupGpuShortRangeWork(gpu_nbv, gpuBonded, iLocality);
360             }
361         }
362
363         //! Returns true if there is GPU short-range work for the given atom locality.
364         bool haveGpuShortRangeWork(const Nbnxm::AtomLocality aLocality)
365         {
366             return ((useGpu() && !emulateGpu()) &&
367                     Nbnxm::haveGpuShortRangeWork(gpu_nbv, aLocality));
368         }
369
370         // TODO: Make all data members private
371     public:
372         //! All data related to the pair lists
373         std::unique_ptr<PairlistSets>     pairlistSets_;
374         //! Working data for constructing the pairlists
375         std::unique_ptr<PairSearch>       pairSearch_;
376         //! Atom data
377         std::unique_ptr<nbnxn_atomdata_t> nbat;
378     private:
379         //! The non-bonded setup, also affects the pairlist construction kernel
380         Nbnxm::KernelSetup                kernelSetup_;
381     public:
382         //! GPU Nbnxm data, only used with a physical GPU (TODO: use unique_ptr)
383         gmx_nbnxn_gpu_t                  *gpu_nbv;
384 };
385
386 namespace Nbnxm
387 {
388
389 /*! \brief Creates an Nbnxm object */
390 std::unique_ptr<nonbonded_verlet_t>
391 init_nb_verlet(const gmx::MDLogger     &mdlog,
392                gmx_bool                 bFEP_NonBonded,
393                const t_inputrec        *ir,
394                const t_forcerec        *fr,
395                const t_commrec         *cr,
396                const gmx_hw_info_t     &hardwareInfo,
397                const gmx_device_info_t *deviceInfo,
398                const gmx_mtop_t        *mtop,
399                matrix                   box);
400
401 } // namespace Nbnxm
402
403 /*! \brief Put the atoms on the pair search grid.
404  *
405  * Only atoms atomStart to atomEnd in x are put on the grid.
406  * The atom_density is used to determine the grid size.
407  * When atomDensity<=0, the density is determined from atomEnd-atomStart and the corners.
408  * With domain decomposition part of the n particles might have migrated,
409  * but have not been removed yet. This count is given by nmoved.
410  * When move[i] < 0 particle i has migrated and will not be put on the grid.
411  * Without domain decomposition move will be NULL.
412  */
413 void nbnxn_put_on_grid(nonbonded_verlet_t             *nb_verlet,
414                        const matrix                    box,
415                        int                             ddZone,
416                        const rvec                      lowerCorner,
417                        const rvec                      upperCorner,
418                        const gmx::UpdateGroupsCog     *updateGroupsCog,
419                        int                             atomStart,
420                        int                             atomEnd,
421                        real                            atomDensity,
422                        gmx::ArrayRef<const int>        atomInfo,
423                        gmx::ArrayRef<const gmx::RVec>  x,
424                        int                             numAtomsMoved,
425                        const int                      *move);
426
427 /*! \brief As nbnxn_put_on_grid, but for the non-local atoms
428  *
429  * with domain decomposition. Should be called after calling
430  * nbnxn_search_put_on_grid for the local atoms / home zone.
431  */
432 void nbnxn_put_on_grid_nonlocal(nonbonded_verlet_t              *nb_verlet,
433                                 const struct gmx_domdec_zones_t *zones,
434                                 gmx::ArrayRef<const int>         atomInfo,
435                                 gmx::ArrayRef<const gmx::RVec>   x);
436
437 #endif // GMX_NBNXN_NBNXM_H