1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
3 * This file is part of GROMACS.
6 * Written by the Gromacs development team under coordination of
7 * David van der Spoel, Berk Hess, and Erik Lindahl.
9 * This library 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
12 * of the License, or (at your option) any later version.
14 * To help us fund GROMACS development, we humbly ask that you cite
15 * the research papers on the package. Check out http://www.gromacs.org
18 * Gnomes, ROck Monsters And Chili Sauce
28 #include "nb_kernel.h"
33 /* Static data structures to find kernels */
34 static nb_kernel_info_t * kernel_list = NULL;
35 static unsigned int kernel_list_size = 0;
36 static int * kernel_list_hash = NULL;
37 static unsigned int kernel_list_hash_size = 0;
40 nb_kernel_hash_func(const char * arch,
42 const char * elec_mod,
51 hash = gmx_string_hash_func(arch,gmx_string_hash_init);
52 hash = gmx_string_hash_func(elec,hash);
53 hash = gmx_string_hash_func(elec_mod,hash);
54 hash = gmx_string_hash_func(vdw,hash);
55 hash = gmx_string_hash_func(vdw_mod,hash);
56 hash = gmx_string_hash_func(geom,hash);
57 hash = gmx_string_hash_func(other,hash);
58 hash = gmx_string_hash_func(vf,hash);
64 nb_kernel_list_add_kernels(nb_kernel_info_t * new_kernel_list,
65 int new_kernel_list_size)
67 srenew(kernel_list,kernel_list_size+new_kernel_list_size);
68 memcpy(kernel_list+kernel_list_size,new_kernel_list,new_kernel_list_size*sizeof(nb_kernel_info_t));
69 kernel_list_size += new_kernel_list_size;
74 nb_kernel_list_hash_init(void)
79 kernel_list_hash_size = kernel_list_size*5;
80 snew(kernel_list_hash,kernel_list_hash_size);
82 for(i=0;i<kernel_list_hash_size;i++)
84 kernel_list_hash[i] = -1;
86 for(i=0;i<kernel_list_size;i++)
88 index=nb_kernel_hash_func(kernel_list[i].architecture,
89 kernel_list[i].electrostatics,
90 kernel_list[i].electrostatics_modifier,
92 kernel_list[i].vdw_modifier,
93 kernel_list[i].geometry,
95 kernel_list[i].vf) % kernel_list_hash_size;
97 /* Check for collisions and advance if necessary */
98 while( kernel_list_hash[index] != -1 )
100 index = (index+1) % kernel_list_hash_size;
103 kernel_list_hash[index] = i;
109 nb_kernel_list_hash_destroy()
111 sfree(kernel_list_hash);
112 kernel_list_hash = NULL;
113 kernel_list_hash_size = 0;
118 nb_kernel_list_findkernel(FILE * log,
120 const char * electrostatics,
121 const char * electrostatics_modifier,
123 const char * vdw_modifier,
124 const char * geometry,
130 nb_kernel_info_t * kernelinfo_ptr;
132 if(kernel_list_hash_size==0)
137 index=nb_kernel_hash_func(arch,
139 electrostatics_modifier,
144 vf) % kernel_list_hash_size;
146 kernelinfo_ptr = NULL;
147 while( (i=kernel_list_hash[index]) != -1)
149 if(!gmx_strcasecmp_min(kernel_list[i].architecture,arch) &&
150 !gmx_strcasecmp_min(kernel_list[i].electrostatics,electrostatics) &&
151 !gmx_strcasecmp_min(kernel_list[i].electrostatics_modifier,electrostatics_modifier) &&
152 !gmx_strcasecmp_min(kernel_list[i].vdw,vdw) &&
153 !gmx_strcasecmp_min(kernel_list[i].vdw_modifier,vdw_modifier) &&
154 !gmx_strcasecmp_min(kernel_list[i].geometry,geometry) &&
155 !gmx_strcasecmp_min(kernel_list[i].other,other) &&
156 !gmx_strcasecmp_min(kernel_list[i].vf,vf))
158 kernelinfo_ptr = kernel_list+i;
161 index = (index+1) % kernel_list_hash_size;
164 if(log && kernelinfo_ptr!=NULL)
167 "NB kernel %s() with architecture '%s' used for neighborlist with\n"
168 " Elec: '%s', Modifier: '%s'\n"
169 " Vdw: '%s', Modifier: '%s'\n"
170 " Geom: '%s', Other: '%s', Calc: '%s'\n\n",
171 kernelinfo_ptr->kernelname,arch,electrostatics,electrostatics_modifier,
172 vdw,vdw_modifier,geometry,other,vf);
175 /* If we did not find any kernel the pointer will still be NULL */
176 return (kernelinfo_ptr !=NULL) ? kernelinfo_ptr->kernelptr : NULL;