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 * Copyright (c) 2011,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 /* This file is completely threadsafe - keep it that way! */
42 /* #define DEBUG_NNB */
43 #include "gpp_nextnb.h"
46 #include "gromacs/utility/fatalerror.h"
47 #include "gromacs/utility/smalloc.h"
54 bond_sort (const void *a, const void *b)
63 return (sa->aj-sb->aj);
67 return (sa->ai-sb->ai);
72 compare_int (const void * a, const void * b)
74 return ( *(int*)a - *(int*)b );
79 #define prints(str, n, s) __prints(str, n, s)
80 static void __prints(char *str, int n, sortable *s)
86 fprintf(debug, "%s\n", str);
87 fprintf(debug, "Sortables \n");
88 for (i = 0; (i < n); i++)
90 fprintf(debug, "%d\t%d\n", s[i].ai, s[i].aj);
97 #define prints(str, n, s)
100 void init_nnb(t_nextnb *nnb, int nr, int nrex)
109 snew(nnb->nrexcl, nr);
110 for (i = 0; (i < nr); i++)
112 snew(nnb->a[i], nrex+1);
113 snew(nnb->nrexcl[i], nrex+1);
117 static void add_nnb (t_nextnb *nnb, int nre, int i, int j)
119 srenew(nnb->a[i][nre], nnb->nrexcl[i][nre]+1);
120 nnb->a[i][nre][nnb->nrexcl[i][nre]] = j;
121 nnb->nrexcl[i][nre]++;
124 void done_nnb (t_nextnb *nnb)
128 for (i = 0; (i < nnb->nr); i++)
130 for (nre = 0; (nre <= nnb->nrex); nre++)
132 if (nnb->nrexcl[i][nre] > 0)
134 sfree (nnb->a[i][nre]);
137 sfree (nnb->nrexcl[i]);
147 void __print_nnb(t_nextnb *nnb, char *s)
153 fprintf(debug, "%s\n", s);
154 fprintf(debug, "nnb->nr: %d\n", nnb->nr);
155 fprintf(debug, "nnb->nrex: %d\n", nnb->nrex);
156 for (i = 0; (i < nnb->nr); i++)
158 for (j = 0; (j <= nnb->nrex); j++)
160 fprintf(debug, "nrexcl[%d][%d]: %d, excl: ", i, j, nnb->nrexcl[i][j]);
161 for (k = 0; (k < nnb->nrexcl[i][j]); k++)
163 fprintf(debug, "%d, ", nnb->a[i][j][k]);
165 fprintf(debug, "\n");
172 static void nnb2excl(t_nextnb *nnb, t_blocka *excl)
175 int nre, nrx, nrs, nr_of_sortables;
178 srenew(excl->index, nnb->nr+1);
180 for (i = 0; (i < nnb->nr); i++)
182 /* calculate the total number of exclusions for atom i */
184 for (nre = 0; (nre <= nnb->nrex); nre++)
186 nr_of_sortables += nnb->nrexcl[i][nre];
191 fprintf(debug, "nr_of_sortables: %d\n", nr_of_sortables);
193 /* make space for sortable array */
194 snew(s, nr_of_sortables);
196 /* fill the sortable array and sort it */
198 for (nre = 0; (nre <= nnb->nrex); nre++)
200 for (nrx = 0; (nrx < nnb->nrexcl[i][nre]); nrx++)
203 s[nrs].aj = nnb->a[i][nre][nrx];
207 if (nrs != nr_of_sortables)
209 gmx_incons("Generating exclusions");
211 prints("nnb2excl before qsort", nr_of_sortables, s);
212 if (nr_of_sortables > 1)
214 qsort ((void *)s, nr_of_sortables, (size_t)sizeof(s[0]), bond_sort);
215 prints("nnb2excl after qsort", nr_of_sortables, s);
218 /* remove duplicate entries from the list */
220 if (nr_of_sortables > 0)
222 for (j = 1; (j < nr_of_sortables); j++)
224 if ((s[j].ai != s[j-1].ai) || (s[j].aj != s[j-1].aj))
226 s[j_index++] = s[j-1];
229 s[j_index++] = s[j-1];
231 nr_of_sortables = j_index;
232 prints("after rm-double", j_index, s);
234 /* make space for arrays */
235 srenew(excl->a, excl->nra+nr_of_sortables);
237 /* put the sorted exclusions in the target list */
238 for (nrs = 0; (nrs < nr_of_sortables); nrs++)
240 excl->a[excl->nra+nrs] = s[nrs].aj;
242 excl->nra += nr_of_sortables;
243 excl->index[i+1] = excl->nra;
245 /* cleanup temporary space */
250 print_blocka(debug, "Exclusions", "Atom", "Excluded", excl);
254 static void do_gen(int nrbonds, /* total number of bonds in s */
255 sortable *s, /* bidirectional list of bonds */
256 t_nextnb *nnb) /* the tmp storage for excl */
257 /* Assume excl is initalised and s[] contains all bonds bidirectional */
262 for (i = 0; (i < nnb->nr); i++)
264 add_nnb(nnb, 0, i, i);
266 print_nnb(nnb, "After exclude self");
268 /* exclude all the bonded atoms */
271 for (i = 0; (i < nrbonds); i++)
273 add_nnb(nnb, 1, s[i].ai, s[i].aj);
276 print_nnb(nnb, "After exclude bonds");
278 /* for the nr of exclusions per atom */
279 for (n = 1; (n < nnb->nrex); n++)
281 /* now for all atoms */
282 for (i = 0; (i < nnb->nr); i++)
284 /* for all directly bonded atoms of atom i */
285 for (j = 0; (j < nnb->nrexcl[i][1]); j++)
288 /* store the 1st neighbour in nb */
289 nb = nnb->a[i][1][j];
291 /* store all atoms in nb's n-th list into i's n+1-th list */
292 for (k = 0; (k < nnb->nrexcl[nb][n]); k++)
294 if (i != nnb->a[nb][n][k])
296 add_nnb(nnb, n+1, i, nnb->a[nb][n][k]);
302 print_nnb(nnb, "After exclude rest");
306 static void add_b(t_params *bonds, int *nrf, sortable *s)
311 for (i = 0; (i < bonds->nr); i++)
313 ai = bonds->param[i].AI;
314 aj = bonds->param[i].AJ;
315 if ((ai < 0) || (aj < 0))
317 gmx_fatal(FARGS, "Impossible atom numbers in bond %d: ai=%d, aj=%d",
320 /* Add every bond twice */
328 void gen_nnb(t_nextnb *nnb, t_params plist[])
334 for (i = 0; (i < F_NRE); i++)
338 /* we need every bond twice (bidirectional) */
339 nrbonds += 2*plist[i].nr;
346 for (i = 0; (i < F_NRE); i++)
350 add_b(&plist[i], &nrf, s);
354 /* now sort the bonds */
355 prints("gen_excl before qsort", nrbonds, s);
358 qsort((void *) s, nrbonds, (size_t)sizeof(sortable), bond_sort);
359 prints("gen_excl after qsort", nrbonds, s);
362 do_gen(nrbonds, s, nnb);
367 sort_and_purge_nnb(t_nextnb *nnb)
369 int i, j, k, m, n, cnt, found, prev, idx;
371 for (i = 0; (i < nnb->nr); i++)
373 for (n = 0; (n <= nnb->nrex); n++)
375 /* Sort atoms in this list */
376 qsort(nnb->a[i][n], nnb->nrexcl[i][n], sizeof(int), compare_int);
380 for (j = 0; j < nnb->nrexcl[i][n]; j++)
382 idx = nnb->a[i][n][j];
385 for (m = 0; m < n && !found; m++)
387 for (k = 0; k < nnb->nrexcl[i][m] && !found; k++)
389 found = (idx == nnb->a[i][m][k]);
393 if (!found && nnb->a[i][n][j] != prev)
395 nnb->a[i][n][cnt] = nnb->a[i][n][j];
396 prev = nnb->a[i][n][cnt];
400 nnb->nrexcl[i][n] = cnt;
406 void generate_excl (int nrexcl, int nratoms, t_params plist[], t_nextnb *nnb, t_blocka *excl)
412 gmx_fatal(FARGS, "Can't have %d exclusions...", nrexcl);
414 init_nnb(nnb, nratoms, nrexcl);
417 sort_and_purge_nnb(nnb);
418 nnb2excl (nnb, excl);