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