Tests of restrained listed potentials.
[alexxy/gromacs.git] / src / gromacs / mdlib / nb_verlet.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2017,2018, 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_nb_verlet 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 NxN algorithms ("nbnxn_" 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 non-bonded Verlet module
86  * that implements the NxN cluster non-bonded algorithm to efficiently compute
87  * 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_nb_verlet
95  */
96
97
98 #ifndef NB_VERLET_H
99 #define NB_VERLET_H
100
101 #include <memory>
102
103 #include "gromacs/mdlib/nbnxn_gpu_types.h"
104 #include "gromacs/mdlib/nbnxn_pairlist.h"
105
106 //! Help pass GPU-emulation parameters with type safety.
107 enum class EmulateGpuNonbonded : bool
108 {
109     //! Do not emulate GPUs.
110     No,
111     //! Do emulate GPUs.
112     Yes
113 };
114
115
116 /*! \brief Nonbonded NxN kernel types: plain C, CPU SIMD, GPU, GPU emulation */
117 typedef enum
118 {
119     nbnxnkNotSet = 0,
120     nbnxnk4x4_PlainC,
121     nbnxnk4xN_SIMD_4xN,
122     nbnxnk4xN_SIMD_2xNN,
123     nbnxnk8x8x8_GPU,
124     nbnxnk8x8x8_PlainC,
125     nbnxnkNR
126 } nbnxn_kernel_type;
127
128 /*! \brief Return a string identifying the kernel type.
129  *
130  * \param [in] kernel_type   nonbonded kernel types, takes values from the nbnxn_kernel_type enum
131  * \returns                  a string identifying the kernel corresponding to the type passed as argument
132  */
133 const char *lookup_nbnxn_kernel_name(int kernel_type);
134
135 /*! \brief Ewald exclusion types */
136 enum {
137     ewaldexclTable, ewaldexclAnalytical
138 };
139
140 /*! \brief Atom locality indicator: local, non-local, all.
141  *
142  * Used for calls to:
143  * gridding, pair-search, force calculation, x/f buffer operations
144  * */
145 enum {
146     eatLocal = 0, eatNonlocal = 1, eatAll
147 };
148
149 /*! \brief Tests for local atom range */
150 #define LOCAL_A(x)               ((x) == eatLocal)
151 /*! \brief Tests for non-local atom range */
152 #define NONLOCAL_A(x)            ((x) == eatNonlocal)
153 /*! \brief Tests for either local or non-local atom range */
154 #define LOCAL_OR_NONLOCAL_A(x)   (LOCAL_A(x) || NONLOCAL_A(x))
155
156 /*! \brief Interaction locality indicator
157  *
158  * Used in pair-list search/calculations in the following manner:
159  *  - local interactions require local atom data and affect local output only;
160  *  - non-local interactions require both local and non-local atom data and
161  *    affect both local- and non-local output.
162  */
163 enum {
164     eintLocal = 0, eintNonlocal = 1
165 };
166
167 /*! \brief Tests for local interaction indicator */
168 #define LOCAL_I(x)               ((x) == eintLocal)
169 /*! \brief Tests for non-local interaction indicator */
170 #define NONLOCAL_I(x)            ((x) == eintNonlocal)
171
172 /*! \brief Flag to tell the nonbonded kernels whether to clear the force output buffers */
173 enum {
174     enbvClearFNo, enbvClearFYes
175 };
176
177 /*! \libinternal
178  *  \brief Non-bonded interaction group data structure. */
179 typedef struct nonbonded_verlet_group_t {
180     nbnxn_pairlist_set_t  nbl_lists;   /**< pair list(s)                       */
181     int                   kernel_type; /**< non-bonded kernel - see enum above */
182     int                   ewald_excl;  /**< Ewald exclusion - see enum above   */
183 } nonbonded_verlet_group_t;
184
185 /*! \libinternal
186  *  \brief Top-level non-bonded data structure for the Verlet-type cut-off scheme. */
187 typedef struct nonbonded_verlet_t {
188     std::unique_ptr<NbnxnListParameters> listParams;      /**< Parameters for the search and list pruning setup */
189     std::unique_ptr<nbnxn_search>        nbs;             /**< n vs n atom pair searching data       */
190     int                                  ngrp;            /**< number of interaction groups          */
191     nonbonded_verlet_group_t             grp[2];          /**< local and non-local interaction group */
192     nbnxn_atomdata_t                    *nbat;            /**< atom data                             */
193
194     gmx_bool                             bUseGPU;         /**< TRUE when non-bonded interactions are computed on a physical GPU */
195     EmulateGpuNonbonded                  emulateGpu;      /**< true when non-bonded interactions are computed on the CPU using GPU-style pair lists */
196     gmx_nbnxn_gpu_t                     *gpu_nbv;         /**< pointer to GPU nb verlet data     */
197     int                                  min_ci_balanced; /**< pair list balancing parameter
198                                                                used for the 8x8x8 GPU kernels    */
199 } nonbonded_verlet_t;
200
201 /*! \brief Getter for bUseGPU */
202 gmx_bool
203 usingGpu(nonbonded_verlet_t *nbv);
204
205 #endif /* NB_VERLET_H */