Tests of restrained listed potentials.
[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/nbnxm/pairlist.h"
105 #include "gromacs/nbnxm/pairlistset.h"
106 #include "gromacs/utility/arrayref.h"
107 #include "gromacs/utility/real.h"
108
109 // TODO: Remove these includes
110 #include "kerneldispatch.h"
111 #include "nbnxm_gpu.h"
112
113 struct gmx_device_info_t;
114 struct gmx_domdec_zones_t;
115 struct gmx_hw_info_t;
116 struct gmx_mtop_t;
117 struct t_commrec;
118 struct t_forcerec;
119 struct t_inputrec;
120
121 namespace gmx
122 {
123 class MDLogger;
124 class UpdateGroupsCog;
125 }
126
127 //! Help pass GPU-emulation parameters with type safety.
128 enum class EmulateGpuNonbonded : bool
129 {
130     //! Do not emulate GPUs.
131     No,
132     //! Do emulate GPUs.
133     Yes
134 };
135
136
137 /*! \brief Nonbonded NxN kernel types: plain C, CPU SIMD, GPU, GPU emulation */
138 typedef enum
139 {
140     nbnxnkNotSet = 0,
141     nbnxnk4x4_PlainC,
142     nbnxnk4xN_SIMD_4xN,
143     nbnxnk4xN_SIMD_2xNN,
144     nbnxnk8x8x8_GPU,
145     nbnxnk8x8x8_PlainC,
146     nbnxnkNR
147 } nbnxn_kernel_type;
148
149 /*! \brief Return a string identifying the kernel type.
150  *
151  * \param [in] kernel_type   nonbonded kernel types, takes values from the nbnxn_kernel_type enum
152  * \returns                  a string identifying the kernel corresponding to the type passed as argument
153  */
154 const char *lookup_nbnxn_kernel_name(int kernel_type);
155
156 /*! \brief Ewald exclusion types */
157 enum {
158     ewaldexclTable, ewaldexclAnalytical
159 };
160
161 /*! \brief Atom locality indicator: local, non-local, all.
162  *
163  * Used for calls to:
164  * gridding, pair-search, force calculation, x/f buffer operations
165  * */
166 enum {
167     eatLocal = 0, eatNonlocal = 1, eatAll
168 };
169
170 /*! \brief Tests for local atom range */
171 #define LOCAL_A(x)               ((x) == eatLocal)
172 /*! \brief Tests for non-local atom range */
173 #define NONLOCAL_A(x)            ((x) == eatNonlocal)
174 /*! \brief Tests for either local or non-local atom range */
175 #define LOCAL_OR_NONLOCAL_A(x)   (LOCAL_A(x) || NONLOCAL_A(x))
176
177 /*! \brief Interaction locality indicator
178  *
179  * Used in pair-list search/calculations in the following manner:
180  *  - local interactions require local atom data and affect local output only;
181  *  - non-local interactions require both local and non-local atom data and
182  *    affect both local- and non-local output.
183  */
184 enum {
185     eintLocal = 0, eintNonlocal = 1
186 };
187
188 /*! \brief Tests for local interaction indicator */
189 #define LOCAL_I(x)               ((x) == eintLocal)
190 /*! \brief Tests for non-local interaction indicator */
191 #define NONLOCAL_I(x)            ((x) == eintNonlocal)
192
193 /*! \brief Flag to tell the nonbonded kernels whether to clear the force output buffers */
194 enum {
195     enbvClearFNo, enbvClearFYes
196 };
197
198 /*! \libinternal
199  *  \brief Non-bonded interaction group data structure. */
200 typedef struct nonbonded_verlet_group_t {
201     nbnxn_pairlist_set_t  nbl_lists;   /**< pair list(s)                       */
202     int                   kernel_type; /**< non-bonded kernel - see enum above */
203     int                   ewald_excl;  /**< Ewald exclusion - see enum above   */
204 } nonbonded_verlet_group_t;
205
206 /*! \libinternal
207  *  \brief Top-level non-bonded data structure for the Verlet-type cut-off scheme. */
208 typedef struct nonbonded_verlet_t {
209     std::unique_ptr<NbnxnListParameters> listParams;      /**< Parameters for the search and list pruning setup */
210     std::unique_ptr<nbnxn_search>        nbs;             /**< n vs n atom pair searching data       */
211     int                                  ngrp;            /**< number of interaction groups          */
212     nonbonded_verlet_group_t             grp[2];          /**< local and non-local interaction group */
213     nbnxn_atomdata_t                    *nbat;            /**< atom data                             */
214
215     gmx_bool                             bUseGPU;         /**< TRUE when non-bonded interactions are computed on a physical GPU */
216     EmulateGpuNonbonded                  emulateGpu;      /**< true when non-bonded interactions are computed on the CPU using GPU-style pair lists */
217     gmx_nbnxn_gpu_t                     *gpu_nbv;         /**< pointer to GPU nb verlet data     */
218     int                                  min_ci_balanced; /**< pair list balancing parameter
219                                                                used for the 8x8x8 GPU kernels    */
220 } nonbonded_verlet_t;
221
222 /*! \brief Initializes the nbnxn module */
223 void init_nb_verlet(const gmx::MDLogger     &mdlog,
224                     nonbonded_verlet_t     **nb_verlet,
225                     gmx_bool                 bFEP_NonBonded,
226                     const t_inputrec        *ir,
227                     const t_forcerec        *fr,
228                     const t_commrec         *cr,
229                     const gmx_hw_info_t     &hardwareInfo,
230                     const gmx_device_info_t *deviceInfo,
231                     const gmx_mtop_t        *mtop,
232                     matrix                   box);
233
234 /*! \brief Put the atoms on the pair search grid.
235  *
236  * Only atoms atomStart to atomEnd in x are put on the grid.
237  * The atom_density is used to determine the grid size.
238  * When atomDensity<=0, the density is determined from atomEnd-atomStart and the corners.
239  * With domain decomposition part of the n particles might have migrated,
240  * but have not been removed yet. This count is given by nmoved.
241  * When move[i] < 0 particle i has migrated and will not be put on the grid.
242  * Without domain decomposition move will be NULL.
243  */
244 void nbnxn_put_on_grid(nbnxn_search_t                  nbs,
245                        int                             ePBC,
246                        const matrix                    box,
247                        int                             ddZone,
248                        const rvec                      lowerCorner,
249                        const rvec                      upperCorner,
250                        const gmx::UpdateGroupsCog     *updateGroupsCog,
251                        int                             atomStart,
252                        int                             atomEnd,
253                        real                            atomDensity,
254                        const int                      *atinfo,
255                        gmx::ArrayRef<const gmx::RVec>  x,
256                        int                             numAtomsMoved,
257                        const int                      *move,
258                        int                             nb_kernel_type,
259                        nbnxn_atomdata_t               *nbat);
260
261 /*! \brief As nbnxn_put_on_grid, but for the non-local atoms
262  *
263  * with domain decomposition. Should be called after calling
264  * nbnxn_search_put_on_grid for the local atoms / home zone.
265  */
266 void nbnxn_put_on_grid_nonlocal(nbnxn_search_t                   nbs,
267                                 const struct gmx_domdec_zones_t *zones,
268                                 const int                       *atinfo,
269                                 gmx::ArrayRef<const gmx::RVec>   x,
270                                 int                              nb_kernel_type,
271                                 nbnxn_atomdata_t                *nbat);
272
273 /*! \brief Returns the number of x and y cells in the local grid */
274 void nbnxn_get_ncells(nbnxn_search_t nbs, int *ncx, int *ncy);
275
276 /*! \brief Returns the order indices of the atoms on the pairlist search grid */
277 gmx::ArrayRef<const int> nbnxn_get_atomorder(const nbnxn_search* nbs);
278
279 /*! \brief Renumbers the atom indices on the grid to consecutive order */
280 void nbnxn_set_atomorder(nbnxn_search_t nbs);
281
282 /*! \brief Returns the index position of the atoms on the pairlist search grid */
283 gmx::ArrayRef<const int> nbnxn_get_gridindices(const nbnxn_search* nbs);
284
285 #endif // GMX_NBNXN_NBNXN_H