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