158dc5eef82d694dc9ad5b4f44a00ac6d57af85a
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,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 #include "gmxpre.h"
36
37 #include "config.h"
38
39
40 #include <stdio.h>
41 #include <string.h>
42
43 #include "nb_kernel.h"
44 #include "gromacs/utility/smalloc.h"
45 #include "gromacs/utility/cstringutil.h"
46 #include "gromacs/utility/fatalerror.h"
47
48
49 /* Static data structures to find kernels */
50 static nb_kernel_info_t *   kernel_list           = NULL;
51 static unsigned int         kernel_list_size      = 0;
52 static int *                kernel_list_hash      = NULL;
53 static unsigned int         kernel_list_hash_size = 0;
54
55 static unsigned int
56 nb_kernel_hash_func(const char *   arch,
57                     const char *   elec,
58                     const char *   elec_mod,
59                     const char *   vdw,
60                     const char *   vdw_mod,
61                     const char *   geom,
62                     const char *   other,
63                     const char *   vf)
64 {
65     unsigned int hash;
66
67     hash = gmx_string_hash_func(arch, gmx_string_hash_init);
68     hash = gmx_string_hash_func(elec, hash);
69     hash = gmx_string_hash_func(elec_mod, hash);
70     hash = gmx_string_hash_func(vdw, hash);
71     hash = gmx_string_hash_func(vdw_mod, hash);
72     hash = gmx_string_hash_func(geom, hash);
73     hash = gmx_string_hash_func(other, hash);
74     hash = gmx_string_hash_func(vf, hash);
75
76     return hash;
77 }
78
79 void
80 nb_kernel_list_add_kernels(nb_kernel_info_t *   new_kernel_list,
81                            int                  new_kernel_list_size)
82 {
83     srenew(kernel_list, kernel_list_size+new_kernel_list_size);
84     memcpy(kernel_list+kernel_list_size, new_kernel_list, new_kernel_list_size*sizeof(nb_kernel_info_t));
85     kernel_list_size += new_kernel_list_size;
86 }
87
88
89 int
90 nb_kernel_list_hash_init(void)
91 {
92     unsigned int            i;
93     unsigned int            index;
94
95     kernel_list_hash_size   = kernel_list_size*5;
96     snew(kernel_list_hash, kernel_list_hash_size);
97
98     for (i = 0; i < kernel_list_hash_size; i++)
99     {
100         kernel_list_hash[i] = -1;
101     }
102     for (i = 0; i < kernel_list_size; i++)
103     {
104         index = nb_kernel_hash_func(kernel_list[i].architecture,
105                                     kernel_list[i].electrostatics,
106                                     kernel_list[i].electrostatics_modifier,
107                                     kernel_list[i].vdw,
108                                     kernel_list[i].vdw_modifier,
109                                     kernel_list[i].geometry,
110                                     kernel_list[i].other,
111                                     kernel_list[i].vf) % kernel_list_hash_size;
112
113         /* Check for collisions and advance if necessary */
114         while (kernel_list_hash[index] != -1)
115         {
116             index = (index+1) % kernel_list_hash_size;
117         }
118
119         kernel_list_hash[index] = i;
120     }
121     return 0;
122 }
123
124 void
125 nb_kernel_list_hash_destroy()
126 {
127     sfree(kernel_list_hash);
128     kernel_list_hash      = NULL;
129     kernel_list_hash_size = 0;
130 }
131
132
133 nb_kernel_t *
134 nb_kernel_list_findkernel(FILE gmx_unused *   log,
135                           const char *        arch,
136                           const char *        electrostatics,
137                           const char *        electrostatics_modifier,
138                           const char *        vdw,
139                           const char *        vdw_modifier,
140                           const char *        geometry,
141                           const char *        other,
142                           const char *        vf)
143 {
144     int                 i;
145     unsigned int        index;
146     nb_kernel_info_t *  kernelinfo_ptr;
147
148     if (kernel_list_hash_size == 0)
149     {
150         return NULL;
151     }
152
153     index = nb_kernel_hash_func(arch,
154                                 electrostatics,
155                                 electrostatics_modifier,
156                                 vdw,
157                                 vdw_modifier,
158                                 geometry,
159                                 other,
160                                 vf) % kernel_list_hash_size;
161
162     kernelinfo_ptr = NULL;
163     while ( (i = kernel_list_hash[index]) != -1)
164     {
165         if (!gmx_strcasecmp_min(kernel_list[i].architecture, arch) &&
166             !gmx_strcasecmp_min(kernel_list[i].electrostatics, electrostatics) &&
167             !gmx_strcasecmp_min(kernel_list[i].electrostatics_modifier, electrostatics_modifier) &&
168             !gmx_strcasecmp_min(kernel_list[i].vdw, vdw) &&
169             !gmx_strcasecmp_min(kernel_list[i].vdw_modifier, vdw_modifier) &&
170             !gmx_strcasecmp_min(kernel_list[i].geometry, geometry) &&
171             !gmx_strcasecmp_min(kernel_list[i].other, other) &&
172             !gmx_strcasecmp_min(kernel_list[i].vf, vf))
173         {
174             kernelinfo_ptr = kernel_list+i;
175             break;
176         }
177         index = (index+1) % kernel_list_hash_size;
178     }
179
180     if (debug && kernelinfo_ptr != NULL)
181     {
182         fprintf(debug,
183                 "NB kernel %s() with architecture '%s' used for neighborlist with\n"
184                 "    Elec: '%s', Modifier: '%s'\n"
185                 "    Vdw:  '%s', Modifier: '%s'\n"
186                 "    Geom: '%s', Other: '%s', Calc: '%s'\n\n",
187                 kernelinfo_ptr->kernelname, arch, electrostatics, electrostatics_modifier,
188                 vdw, vdw_modifier, geometry, other, vf);
189     }
190
191     /* If we did not find any kernel the pointer will still be NULL */
192     return (kernelinfo_ptr != NULL) ? kernelinfo_ptr->kernelptr : NULL;
193 }