4c84070ddedc6d14684bc464137e0bbd8451853c
[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 enum class PairlistType;
122 class PairSearch;
123 class PairlistSet;
124 struct t_blocka;
125 struct t_commrec;
126 struct t_lambda;
127 struct t_mdatoms;
128 struct t_nrnb;
129 struct t_forcerec;
130 struct t_inputrec;
131
132 namespace gmx
133 {
134 class MDLogger;
135 class UpdateGroupsCog;
136 }
137
138 namespace Nbnxm
139 {
140 enum class KernelType;
141 }
142
143 /*! \libinternal
144  * \brief The setup for generating and pruning the nbnxn pair list.
145  *
146  * Without dynamic pruning rlistOuter=rlistInner.
147  */
148 struct NbnxnListParameters
149 {
150     /*! \brief Constructor producing a struct with dynamic pruning disabled
151      */
152     NbnxnListParameters(Nbnxm::KernelType kernelType,
153                         bool              haveFep,
154                         real              rlist,
155                         bool              haveMultipleDomains);
156
157     PairlistType pairlistType;           //!< The type of cluster-pair list
158     bool         haveFep;                //!< Tells whether we have perturbed interactions
159     real         rlistOuter;             //!< Cut-off of the larger, outer pair-list
160     real         rlistInner;             //!< Cut-off of the smaller, inner pair-list
161     bool         haveMultipleDomains;    //!< True when using DD with multiple domains
162     bool         useDynamicPruning;      //!< Are we using dynamic pair-list pruning
163     int          nstlistPrune;           //!< Pair-list dynamic pruning interval
164     int          numRollingPruningParts; //!< The number parts to divide the pair-list into for rolling pruning, a value of 1 gives no rolling pruning
165     int          lifetime;               //!< Lifetime in steps of the pair-list
166 };
167
168 /*! \brief Resources that can be used to execute non-bonded kernels on */
169 enum class NonbondedResource : int
170 {
171     Cpu,
172     Gpu,
173     EmulateGpu
174 };
175
176 namespace Nbnxm
177 {
178
179 /*! \brief Nonbonded NxN kernel types: plain C, CPU SIMD, GPU, GPU emulation */
180 enum class KernelType : int
181 {
182     NotSet = 0,
183     Cpu4x4_PlainC,
184     Cpu4xN_Simd_4xN,
185     Cpu4xN_Simd_2xNN,
186     Gpu8x8x8,
187     Cpu8x8x8_PlainC,
188     Count
189 };
190
191 /*! \brief Ewald exclusion types */
192 enum class EwaldExclusionType : int
193 {
194     NotSet = 0,
195     Table,
196     Analytical,
197     DecidedByGpuModule
198 };
199
200 /* \brief The non-bonded setup, also affects the pairlist construction kernel */
201 struct KernelSetup
202 {
203     //! The non-bonded type, also affects the pairlist construction kernel
204     KernelType         kernelType = KernelType::NotSet;
205     //! Ewald exclusion computation handling type, currently only used for CPU
206     EwaldExclusionType ewaldExclusionType = EwaldExclusionType::NotSet;
207 };
208
209 /*! \brief Return a string identifying the kernel type.
210  *
211  * \param [in] kernelType   nonbonded kernel type, takes values from the nbnxn_kernel_type enum
212  * \returns                 a string identifying the kernel corresponding to the type passed as argument
213  */
214 const char *lookup_kernel_name(Nbnxm::KernelType kernelType);
215
216 } // namespace Nbnxm
217
218 /*! \brief Flag to tell the nonbonded kernels whether to clear the force output buffers */
219 enum {
220     enbvClearFNo, enbvClearFYes
221 };
222
223 /*! \brief Generates a pair-list for the given locality.
224  *
225  * With perturbed particles, also a group scheme style nbl_fep list is made.
226  */
227 void nbnxn_make_pairlist(nonbonded_verlet_t         *nbv,
228                          Nbnxm::InteractionLocality  iLocality,
229                          PairlistSet                *pairlistSet,
230                          const t_blocka             *excl,
231                          int64_t                     step,
232                          t_nrnb                     *nrnb);
233
234 /*! \brief Prune all pair-lists with given locality (currently CPU only)
235  *
236  * For all pair-lists with given locality, takes the outer list and prunes out
237  * pairs beyond the pairlist inner radius and writes the result to a list that is
238  * to be consumed by the non-bonded kernel.
239  */
240 void NbnxnDispatchPruneKernel(PairlistSet            *pairlistSet,
241                               Nbnxm::KernelType       kernelType,
242                               const nbnxn_atomdata_t *nbat,
243                               const rvec             *shift_vec);
244
245 /*! \libinternal
246  *  \brief Top-level non-bonded data structure for the Verlet-type cut-off scheme. */
247 struct nonbonded_verlet_t
248 {
249     public:
250         class PairlistSets
251         {
252             public:
253                 PairlistSets(const NbnxnListParameters  &listParams,
254                              bool                        haveMultipleDomains,
255                              int                         minimumIlistCountForGpuBalancing);
256
257                 //! Construct the pairlist set for the given locality
258                 void construct(Nbnxm::InteractionLocality  iLocality,
259                                PairSearch                 *pairSearch,
260                                nbnxn_atomdata_t           *nbat,
261                                const t_blocka             *excl,
262                                Nbnxm::KernelType           kernelbType,
263                                int64_t                     step,
264                                t_nrnb                     *nrnb);
265
266                 //! Dispatches the dynamic pruning kernel for the given locality
267                 void dispatchPruneKernel(Nbnxm::InteractionLocality  iLocality,
268                                          const nbnxn_atomdata_t     *nbat,
269                                          const rvec                 *shift_vec,
270                                          Nbnxm::KernelType           kernelbType);
271
272                 //! Returns the pair list parameters
273                 const NbnxnListParameters &params() const
274                 {
275                     return params_;
276                 }
277
278                 //! Returns the number of steps performed with the current pair list
279                 int numStepsWithPairlist(int64_t step) const
280                 {
281                     return step - outerListCreationStep_;
282                 }
283
284                 //! Returns whether step is a dynamic list pruning step, for CPU lists
285                 bool isDynamicPruningStepCpu(int64_t step) const
286                 {
287                     return (params_.useDynamicPruning &&
288                             numStepsWithPairlist(step) % params_.nstlistPrune == 0);
289                 }
290
291                 //! Returns whether step is a dynamic list pruning step, for GPU lists
292                 bool isDynamicPruningStepGpu(int64_t step) const
293                 {
294                     const int age = numStepsWithPairlist(step);
295
296                     return (params_.useDynamicPruning &&
297                             age > 0 &&
298                             age < params_.lifetime &&
299                             (params_.haveMultipleDomains || age % 2 == 0));
300                 }
301
302                 //! Changes the pair-list outer and inner radius
303                 void changeRadii(real rlistOuter,
304                                  real rlistInner)
305                 {
306                     params_.rlistOuter = rlistOuter;
307                     params_.rlistInner = rlistInner;
308                 }
309
310                 //! Returns the pair-list set for the given locality
311                 const PairlistSet &pairlistSet(Nbnxm::InteractionLocality iLocality) const
312                 {
313                     if (iLocality == Nbnxm::InteractionLocality::Local)
314                     {
315                         return *localSet_;
316                     }
317                     else
318                     {
319                         GMX_ASSERT(nonlocalSet_, "Need a non-local set when requesting access");
320                         return *nonlocalSet_;
321                     }
322                 }
323
324             private:
325                 //! Returns the pair-list set for the given locality
326                 PairlistSet &pairlistSet(Nbnxm::InteractionLocality iLocality)
327                 {
328                     if (iLocality == Nbnxm::InteractionLocality::Local)
329                     {
330                         return *localSet_;
331                     }
332                     else
333                     {
334                         GMX_ASSERT(nonlocalSet_, "Need a non-local set when requesting access");
335                         return *nonlocalSet_;
336                     }
337                 }
338
339                 //! Parameters for the search and list pruning setup
340                 NbnxnListParameters                   params_;
341                 //! Pair list balancing parameter for use with GPU
342                 int                                   minimumIlistCountForGpuBalancing_;
343                 //! Local pairlist set
344                 std::unique_ptr<PairlistSet>          localSet_;
345                 //! Non-local pairlist set
346                 std::unique_ptr<PairlistSet>          nonlocalSet_;
347                 //! MD step at with the outer lists in pairlistSets_ were created
348                 int64_t                               outerListCreationStep_;
349         };
350
351         //! Constructs an object from its components
352         nonbonded_verlet_t(std::unique_ptr<PairlistSets>      pairlistSets,
353                            std::unique_ptr<PairSearch>        pairSearch,
354                            std::unique_ptr<nbnxn_atomdata_t>  nbat,
355                            const Nbnxm::KernelSetup          &kernelSetup,
356                            gmx_nbnxn_gpu_t                   *gpu_nbv);
357
358         ~nonbonded_verlet_t();
359
360         //! Returns whether a GPU is use for the non-bonded calculations
361         bool useGpu() const
362         {
363             return kernelSetup_.kernelType == Nbnxm::KernelType::Gpu8x8x8;
364         }
365
366         //! Returns whether a GPU is emulated for the non-bonded calculations
367         bool emulateGpu() const
368         {
369             return kernelSetup_.kernelType == Nbnxm::KernelType::Cpu8x8x8_PlainC;
370         }
371
372         //! Return whether the pairlist is of simple, CPU type
373         bool pairlistIsSimple() const
374         {
375             return !useGpu() && !emulateGpu();
376         }
377
378         //! Initialize the pair list sets, TODO this should be private
379         void initPairlistSets(bool haveMultipleDomains);
380
381         //! Returns the order of the local atoms on the grid
382         gmx::ArrayRef<const int> getLocalAtomOrder() const;
383
384         //! Sets the order of the local atoms to the order grid atom ordering
385         void setLocalAtomOrder();
386
387         //! Returns the index position of the atoms on the search grid
388         gmx::ArrayRef<const int> getGridIndices() const;
389
390         //! Constructs the pairlist for the given locality
391         void constructPairlist(Nbnxm::InteractionLocality  iLocality,
392                                const t_blocka             *excl,
393                                int64_t                     step,
394                                t_nrnb                     *nrnb);
395
396         //! Updates all the atom properties in Nbnxm
397         void setAtomProperties(const t_mdatoms &mdatoms,
398                                const int       &atinfo);
399
400         //! Updates the coordinates in Nbnxm for the given locality
401         void setCoordinates(Nbnxm::AtomLocality             locality,
402                             bool                            fillLocal,
403                             gmx::ArrayRef<const gmx::RVec>  x,
404                             gmx_wallcycle                  *wcycle);
405
406         //! Returns a reference to the pairlist sets
407         const PairlistSets &pairlistSets() const
408         {
409             return *pairlistSets_;
410         }
411
412         //! Dispatches the dynamic pruning kernel for the given locality, for CPU lists
413         void dispatchPruneKernelCpu(Nbnxm::InteractionLocality  iLocality,
414                                     const rvec                 *shift_vec);
415
416         //! Dispatches the dynamic pruning kernel for GPU lists
417         void dispatchPruneKernelGpu(int64_t step)
418         {
419             const bool stepIsEven = (pairlistSets().numStepsWithPairlist(step) % 2 == 0);
420
421             Nbnxm::gpu_launch_kernel_pruneonly(gpu_nbv,
422                                                stepIsEven ? Nbnxm::InteractionLocality::Local : Nbnxm::InteractionLocality::NonLocal,
423                                                pairlistSets().params().numRollingPruningParts);
424         }
425
426         //! \brief Executes the non-bonded kernel of the GPU or launches it on the GPU
427         void dispatchNonbondedKernel(Nbnxm::InteractionLocality  iLocality,
428                                      const interaction_const_t  &ic,
429                                      int                         forceFlags,
430                                      int                         clearF,
431                                      t_forcerec                 *fr,
432                                      gmx_enerdata_t             *enerd,
433                                      t_nrnb                     *nrnb);
434
435         //! Executes the non-bonded free-energy kernel, always runs on the CPU
436         void dispatchFreeEnergyKernel(Nbnxm::InteractionLocality  iLocality,
437                                       t_forcerec                 *fr,
438                                       rvec                        x[],
439                                       rvec                        f[],
440                                       const t_mdatoms            &mdatoms,
441                                       t_lambda                   *fepvals,
442                                       real                       *lambda,
443                                       gmx_enerdata_t             *enerd,
444                                       int                         forceFlags,
445                                       t_nrnb                     *nrnb);
446
447         //! Add the forces stored in nbat to f, zeros the forces in nbat */
448         void atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality  locality,
449                                       rvec                *f,
450                                       gmx_wallcycle       *wcycle);
451
452         //! Return the kernel setup
453         const Nbnxm::KernelSetup &kernelSetup() const
454         {
455             return kernelSetup_;
456         }
457
458         /*! \brief Returns the number of x and y cells in the local grid */
459         void getLocalNumCells(int *numCellsX,
460                               int *numCellsY) const;
461
462         // TODO: Make all data members private
463     public:
464         //! All data related to the pair lists
465         std::unique_ptr<PairlistSets>     pairlistSets_;
466         //! Working data for constructing the pairlists
467         std::unique_ptr<PairSearch>       pairSearch_;
468         //! Atom data
469         std::unique_ptr<nbnxn_atomdata_t> nbat;
470     private:
471         //! The non-bonded setup, also affects the pairlist construction kernel
472         Nbnxm::KernelSetup                kernelSetup_;
473     public:
474         //! GPU Nbnxm data, only used with a physical GPU (TODO: use unique_ptr)
475         gmx_nbnxn_gpu_t                  *gpu_nbv;
476 };
477
478 namespace Nbnxm
479 {
480
481 /*! \brief Creates an Nbnxm object */
482 std::unique_ptr<nonbonded_verlet_t>
483 init_nb_verlet(const gmx::MDLogger     &mdlog,
484                gmx_bool                 bFEP_NonBonded,
485                const t_inputrec        *ir,
486                const t_forcerec        *fr,
487                const t_commrec         *cr,
488                const gmx_hw_info_t     &hardwareInfo,
489                const gmx_device_info_t *deviceInfo,
490                const gmx_mtop_t        *mtop,
491                matrix                   box);
492
493 } // namespace Nbnxm
494
495 /*! \brief Put the atoms on the pair search grid.
496  *
497  * Only atoms atomStart to atomEnd in x are put on the grid.
498  * The atom_density is used to determine the grid size.
499  * When atomDensity<=0, the density is determined from atomEnd-atomStart and the corners.
500  * With domain decomposition part of the n particles might have migrated,
501  * but have not been removed yet. This count is given by nmoved.
502  * When move[i] < 0 particle i has migrated and will not be put on the grid.
503  * Without domain decomposition move will be NULL.
504  */
505 void nbnxn_put_on_grid(nonbonded_verlet_t             *nb_verlet,
506                        const matrix                    box,
507                        int                             ddZone,
508                        const rvec                      lowerCorner,
509                        const rvec                      upperCorner,
510                        const gmx::UpdateGroupsCog     *updateGroupsCog,
511                        int                             atomStart,
512                        int                             atomEnd,
513                        real                            atomDensity,
514                        const int                      *atinfo,
515                        gmx::ArrayRef<const gmx::RVec>  x,
516                        int                             numAtomsMoved,
517                        const int                      *move);
518
519 /*! \brief As nbnxn_put_on_grid, but for the non-local atoms
520  *
521  * with domain decomposition. Should be called after calling
522  * nbnxn_search_put_on_grid for the local atoms / home zone.
523  */
524 void nbnxn_put_on_grid_nonlocal(nonbonded_verlet_t              *nb_verlet,
525                                 const struct gmx_domdec_zones_t *zones,
526                                 const int                       *atinfo,
527                                 gmx::ArrayRef<const gmx::RVec>   x);
528
529 #endif // GMX_NBNXN_NBNXM_H