643ac97af7100b1bd80b21315af40e726619c46c
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,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 #ifndef _nb_kernel_h_
36 #define _nb_kernel_h_
37
38 #include <stdio.h>
39
40 #include "gromacs/gmxlib/nrnb.h"
41 #include "gromacs/math/vectypes.h"
42 #include "gromacs/mdtypes/forcerec.h"
43 #include "gromacs/mdtypes/mdatom.h"
44 #include "gromacs/mdtypes/nblist.h"
45 #include "gromacs/topology/block.h"
46 #include "gromacs/utility/real.h"
47
48 #ifdef __cplusplus
49 extern "C" {
50 #endif
51 #if 0
52 } /* fixes auto-indentation problems */
53 #endif
54
55 /*! \def gmx_inline
56  * \brief
57  * Keyword to use in C code instead of C99 `inline`.
58  *
59  * Some of the C compilers we support do not recognize the C99 keyword
60  * `inline`.  This macro should be used in C code and in shared C/C++ headers
61  * to indicate a function is inlined.
62  * C++ code should use plain `inline`, as that is already in C++98.
63  */
64 #if !defined __cplusplus && defined _MSC_VER
65 #define gmx_inline __inline
66 #else
67 /* C++ or C99 */
68 #define gmx_inline inline
69 #endif
70
71 /* Structure to collect kernel data not available in forcerec or mdatoms structures.
72  * This is only used inside the nonbonded module.
73  */
74 typedef struct
75 {
76     int                flags;
77     const t_blocka    *exclusions;
78     real *             lambda;
79     real *             dvdl;
80
81     /* pointers to tables */
82     t_forcetable *     table_elec;
83     t_forcetable *     table_vdw;
84     t_forcetable *     table_elec_vdw;
85
86     /* potentials */
87     real *             energygrp_elec;
88     real *             energygrp_vdw;
89 }
90 nb_kernel_data_t;
91
92
93 typedef void
94     nb_kernel_t (t_nblist *                nlist,
95                  rvec *                    x,
96                  rvec *                    f,
97                  const struct t_forcerec * fr,
98                  const t_mdatoms *         mdatoms,
99                  nb_kernel_data_t *        kernel_data,
100                  t_nrnb *                  nrnb);
101
102
103 /* Structure with a kernel pointer and settings. This cannot be abstract
104  * since we define the kernel list statically for each architecture in a header,
105  * and use it to set up the kernel hash functions to find kernels.
106  *
107  * The electrostatics/vdw names should be obvious and correspond to the
108  * forms of the interactions calculated in this function, and the interaction
109  * modifiers describe switch/shift and similar alterations. Geometry refers
110  * to whether this kernel calculates interactions between single particles or
111  * waters (groups of 3/4 atoms) for better performance. Finally, the VF string
112  * selects whether the kernel calculates only potential, only force, or both
113  *
114  * The allowed values for kernel interactions are described by the
115  * enumerated types gmx_nbkernel_elec and gmx_nbkernel_vdw (see types/enums.h).
116  * Note that these are deliberately NOT identical to the interactions the
117  * user can set, since some user-specified interactions will be tabulated, and
118  * Lennard-Jones and Buckingham use different kernels while their setting in
119  * the input is decided by nonbonded parameter formats rather than mdp options.
120  *
121  * The interaction modifiers are described by the eintmod enum type, while
122  * the kernel geometry is decided from the neighborlist geometry, which is
123  * described by the enum gmx_nblist_kernel_geometry (again, see types/enums.h).
124  * The
125  *
126  * Note that any particular implementation of kernels might not support all of
127  * these strings. In fact, some might not be supported by any architecture yet.
128  * The whole point of using strings and hashes is that we do not have to define a
129  * unique set of strings in a single place. Thus, as long as you implement a
130  * corresponding kernel, you could in theory provide any string you want.
131  */
132 typedef struct nb_kernel_info
133 {
134     nb_kernel_t *   kernelptr;
135     const char *    kernelname;
136     const char *    architecture;     /* e.g. "C", "SSE", "AVX_256", etc. */
137
138     const char *    electrostatics;
139     const char *    electrostatics_modifier;
140     const char *    vdw;
141     const char *    vdw_modifier;
142     const char *    geometry;
143     const char *    other;  /* Any extra info you want/need to select a kernel */
144     const char *    vf;     /* "PotentialAndForce", "Potential", or "Force" */
145 }
146 nb_kernel_info_t;
147
148
149 void
150 nb_kernel_list_add_kernels(nb_kernel_info_t *   new_kernelinfo,
151                            int                  new_size);
152
153 int
154 nb_kernel_list_hash_init(void);
155
156 /* Return a function pointer to the nonbonded kernel pointer with
157  * settings according to the text strings provided. GROMACS does not guarantee
158  * the existence of accelerated kernels for any combination, so the return value
159  * can be NULL.
160  * In that case, you can try a different/lower-level acceleration, and
161  * eventually you need to prepare to fall back to generic kernels or change
162  * your settings and try again.
163  *
164  * The names of the text strings are obviously meant to reflect settings in
165  * GROMACS, but inside this routine they are merely used as a set of text
166  * strings not defined here. The routine will simply compare the arguments with
167  * the contents of the corresponding strings in the nb_kernel_list_t structure.
168  *
169  * This function does not check whether the kernel in question can run on the
170  * present architecture since that would require a slow cpuid call for every
171  * single invocation.
172  */
173 nb_kernel_t *
174 nb_kernel_list_findkernel(FILE *              log,
175                           const char *        architecture,
176                           const char *        electrostatics,
177                           const char *        electrostatics_modifier,
178                           const char *        vdw,
179                           const char *        vdw_modifier,
180                           const char *        geometry,
181                           const char *        other,
182                           const char *        vf);
183
184
185
186 #ifdef __cplusplus
187 }
188 #endif
189
190 #endif /* _nb_kernel_h_ */