2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 /* This file is completely threadsafe - keep it that way! */
46 /* #define DEBUG_NNB */
47 #include "gpp_nextnb.h"
48 #include "gmx_fatal.h"
56 bond_sort (const void *a, const void *b)
64 return (sa->aj-sb->aj);
66 return (sa->ai-sb->ai);
70 compare_int (const void * a, const void * b)
72 return ( *(int*)a - *(int*)b );
77 #define prints(str, n, s) __prints(str, n, s)
78 static void __prints(char *str, int n, sortable *s)
83 fprintf(debug,"%s\n",str);
84 fprintf(debug,"Sortables \n");
85 for (i=0; (i < n); i++)
86 fprintf(debug,"%d\t%d\n",s[i].ai,s[i].aj);
92 #define prints(str, n, s)
95 void init_nnb(t_nextnb *nnb, int nr, int nrex)
104 snew(nnb->nrexcl,nr);
105 for (i=0; (i < nr); i++) {
106 snew(nnb->a[i],nrex+1);
107 snew(nnb->nrexcl[i],nrex+1);
111 static void add_nnb (t_nextnb *nnb, int nre, int i, int j)
113 srenew(nnb->a[i][nre],nnb->nrexcl[i][nre]+1);
114 nnb->a[i][nre][nnb->nrexcl[i][nre]] = j;
115 nnb->nrexcl[i][nre]++;
118 void done_nnb (t_nextnb *nnb)
122 for (i=0; (i < nnb->nr); i++) {
123 for (nre=0; (nre <= nnb->nrex); nre++)
124 if (nnb->nrexcl[i][nre] > 0)
125 sfree (nnb->a[i][nre]);
126 sfree (nnb->nrexcl[i]);
136 void __print_nnb(t_nextnb *nnb, char *s)
141 fprintf(debug,"%s\n",s);
142 fprintf(debug,"nnb->nr: %d\n",nnb->nr);
143 fprintf(debug,"nnb->nrex: %d\n",nnb->nrex);
144 for(i=0; (i<nnb->nr); i++)
145 for(j=0; (j<=nnb->nrex); j++) {
146 fprintf(debug,"nrexcl[%d][%d]: %d, excl: ",i,j,nnb->nrexcl[i][j]);
147 for(k=0; (k<nnb->nrexcl[i][j]); k++)
148 fprintf(debug,"%d, ",nnb->a[i][j][k]);
155 void nnb2excl (t_nextnb *nnb, t_blocka *excl)
158 int nre,nrx,nrs,nr_of_sortables;
161 srenew(excl->index, nnb->nr+1);
163 for (i=0; (i < nnb->nr); i++) {
164 /* calculate the total number of exclusions for atom i */
166 for (nre=0; (nre<=nnb->nrex); nre++)
167 nr_of_sortables+=nnb->nrexcl[i][nre];
170 fprintf(debug,"nr_of_sortables: %d\n",nr_of_sortables);
171 /* make space for sortable array */
172 snew(s,nr_of_sortables);
174 /* fill the sortable array and sort it */
176 for (nre=0; (nre <= nnb->nrex); nre++)
177 for (nrx=0; (nrx < nnb->nrexcl[i][nre]); nrx++) {
179 s[nrs].aj = nnb->a[i][nre][nrx];
182 if (nrs != nr_of_sortables)
183 gmx_incons("Generating exclusions");
184 prints("nnb2excl before qsort",nr_of_sortables,s);
185 if (nr_of_sortables > 1) {
186 qsort ((void *)s,nr_of_sortables,(size_t)sizeof(s[0]),bond_sort);
187 prints("nnb2excl after qsort",nr_of_sortables,s);
190 /* remove duplicate entries from the list */
192 if (nr_of_sortables > 0) {
193 for (j=1; (j < nr_of_sortables); j++)
194 if ((s[j].ai!=s[j-1].ai) || (s[j].aj!=s[j-1].aj))
198 nr_of_sortables=j_index;
199 prints("after rm-double",j_index,s);
201 /* make space for arrays */
202 srenew(excl->a,excl->nra+nr_of_sortables);
204 /* put the sorted exclusions in the target list */
205 for (nrs=0; (nrs<nr_of_sortables); nrs++)
206 excl->a[excl->nra+nrs] = s[nrs].aj;
207 excl->nra+=nr_of_sortables;
208 excl->index[i+1]=excl->nra;
210 /* cleanup temporary space */
214 print_blocka(debug,"Exclusions","Atom","Excluded",excl);
217 static void do_gen(int nrbonds, /* total number of bonds in s */
218 sortable *s, /* bidirectional list of bonds */
219 t_nextnb *nnb) /* the tmp storage for excl */
220 /* Assume excl is initalised and s[] contains all bonds bidirectional */
225 for(i=0; (i<nnb->nr); i++)
229 print_nnb(nnb,"After exclude self");
231 /* exclude all the bonded atoms */
234 for (i=0; (i < nrbonds); i++)
236 add_nnb(nnb,1,s[i].ai,s[i].aj);
239 print_nnb(nnb,"After exclude bonds");
241 /* for the nr of exclusions per atom */
242 for (n=1; (n < nnb->nrex); n++)
244 /* now for all atoms */
245 for (i=0; (i < nnb->nr); i++)
247 /* for all directly bonded atoms of atom i */
248 for (j=0; (j < nnb->nrexcl[i][1]); j++)
251 /* store the 1st neighbour in nb */
252 nb = nnb->a[i][1][j];
254 /* store all atoms in nb's n-th list into i's n+1-th list */
255 for (k=0; (k < nnb->nrexcl[nb][n]); k++)
257 if (i != nnb->a[nb][n][k])
259 add_nnb(nnb,n+1,i,nnb->a[nb][n][k]);
265 print_nnb(nnb,"After exclude rest");
269 static void add_b(t_params *bonds, int *nrf, sortable *s)
274 for (i=0; (i < bonds->nr); i++) {
275 ai = bonds->param[i].AI;
276 aj = bonds->param[i].AJ;
277 if ((ai < 0) || (aj < 0))
278 gmx_fatal(FARGS,"Impossible atom numbers in bond %d: ai=%d, aj=%d",
280 /* Add every bond twice */
288 void gen_nnb(t_nextnb *nnb,t_params plist[])
294 for(i=0; (i<F_NRE); i++)
296 /* we need every bond twice (bidirectional) */
297 nrbonds += 2*plist[i].nr;
302 for(i=0; (i<F_NRE); i++)
304 add_b(&plist[i],&nrf,s);
306 /* now sort the bonds */
307 prints("gen_excl before qsort",nrbonds,s);
309 qsort((void *) s,nrbonds,(size_t)sizeof(sortable),bond_sort);
310 prints("gen_excl after qsort",nrbonds,s);
313 do_gen(nrbonds,s,nnb);
318 sort_and_purge_nnb(t_nextnb *nnb)
320 int i,j,k,m,n,cnt,found,prev,idx;
322 for (i=0; (i < nnb->nr); i++)
324 for (n=0; (n <= nnb->nrex); n++)
326 /* Sort atoms in this list */
327 qsort(nnb->a[i][n],nnb->nrexcl[i][n],sizeof(int),compare_int);
331 for(j=0;j<nnb->nrexcl[i][n];j++)
333 idx = nnb->a[i][n][j];
336 for(m=0;m<n && !found;m++)
338 for(k=0;k<nnb->nrexcl[i][m] && !found;k++)
340 found = (idx==nnb->a[i][m][k]);
344 if(!found && nnb->a[i][n][j]!=prev)
346 nnb->a[i][n][cnt]=nnb->a[i][n][j];
347 prev = nnb->a[i][n][cnt];
351 nnb->nrexcl[i][n]=cnt;
357 void generate_excl (int nrexcl,int nratoms,t_params plist[],t_nextnb *nnb,t_blocka *excl)
363 gmx_fatal(FARGS,"Can't have %d exclusions...",nrexcl);
365 init_nnb(nnb,nratoms,nrexcl);
368 sort_and_purge_nnb(nnb);