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! */
45 /* #define DEBUG_NNB */
46 #include "gpp_nextnb.h"
47 #include "gmx_fatal.h"
55 bond_sort (const void *a, const void *b)
64 return (sa->aj-sb->aj);
68 return (sa->ai-sb->ai);
73 compare_int (const void * a, const void * b)
75 return ( *(int*)a - *(int*)b );
80 #define prints(str, n, s) __prints(str, n, s)
81 static void __prints(char *str, int n, sortable *s)
87 fprintf(debug, "%s\n", str);
88 fprintf(debug, "Sortables \n");
89 for (i = 0; (i < n); i++)
91 fprintf(debug, "%d\t%d\n", s[i].ai, s[i].aj);
98 #define prints(str, n, s)
101 void init_nnb(t_nextnb *nnb, int nr, int nrex)
110 snew(nnb->nrexcl, nr);
111 for (i = 0; (i < nr); i++)
113 snew(nnb->a[i], nrex+1);
114 snew(nnb->nrexcl[i], nrex+1);
118 static void add_nnb (t_nextnb *nnb, int nre, int i, int j)
120 srenew(nnb->a[i][nre], nnb->nrexcl[i][nre]+1);
121 nnb->a[i][nre][nnb->nrexcl[i][nre]] = j;
122 nnb->nrexcl[i][nre]++;
125 void done_nnb (t_nextnb *nnb)
129 for (i = 0; (i < nnb->nr); i++)
131 for (nre = 0; (nre <= nnb->nrex); nre++)
133 if (nnb->nrexcl[i][nre] > 0)
135 sfree (nnb->a[i][nre]);
138 sfree (nnb->nrexcl[i]);
148 void __print_nnb(t_nextnb *nnb, char *s)
154 fprintf(debug, "%s\n", s);
155 fprintf(debug, "nnb->nr: %d\n", nnb->nr);
156 fprintf(debug, "nnb->nrex: %d\n", nnb->nrex);
157 for (i = 0; (i < nnb->nr); i++)
159 for (j = 0; (j <= nnb->nrex); j++)
161 fprintf(debug, "nrexcl[%d][%d]: %d, excl: ", i, j, nnb->nrexcl[i][j]);
162 for (k = 0; (k < nnb->nrexcl[i][j]); k++)
164 fprintf(debug, "%d, ", nnb->a[i][j][k]);
166 fprintf(debug, "\n");
173 static void nnb2excl(t_nextnb *nnb, t_blocka *excl)
176 int nre, nrx, nrs, nr_of_sortables;
179 srenew(excl->index, nnb->nr+1);
181 for (i = 0; (i < nnb->nr); i++)
183 /* calculate the total number of exclusions for atom i */
185 for (nre = 0; (nre <= nnb->nrex); nre++)
187 nr_of_sortables += nnb->nrexcl[i][nre];
192 fprintf(debug, "nr_of_sortables: %d\n", nr_of_sortables);
194 /* make space for sortable array */
195 snew(s, nr_of_sortables);
197 /* fill the sortable array and sort it */
199 for (nre = 0; (nre <= nnb->nrex); nre++)
201 for (nrx = 0; (nrx < nnb->nrexcl[i][nre]); nrx++)
204 s[nrs].aj = nnb->a[i][nre][nrx];
208 if (nrs != nr_of_sortables)
210 gmx_incons("Generating exclusions");
212 prints("nnb2excl before qsort", nr_of_sortables, s);
213 if (nr_of_sortables > 1)
215 qsort ((void *)s, nr_of_sortables, (size_t)sizeof(s[0]), bond_sort);
216 prints("nnb2excl after qsort", nr_of_sortables, s);
219 /* remove duplicate entries from the list */
221 if (nr_of_sortables > 0)
223 for (j = 1; (j < nr_of_sortables); j++)
225 if ((s[j].ai != s[j-1].ai) || (s[j].aj != s[j-1].aj))
227 s[j_index++] = s[j-1];
230 s[j_index++] = s[j-1];
232 nr_of_sortables = j_index;
233 prints("after rm-double", j_index, s);
235 /* make space for arrays */
236 srenew(excl->a, excl->nra+nr_of_sortables);
238 /* put the sorted exclusions in the target list */
239 for (nrs = 0; (nrs < nr_of_sortables); nrs++)
241 excl->a[excl->nra+nrs] = s[nrs].aj;
243 excl->nra += nr_of_sortables;
244 excl->index[i+1] = excl->nra;
246 /* cleanup temporary space */
251 print_blocka(debug, "Exclusions", "Atom", "Excluded", excl);
255 static void do_gen(int nrbonds, /* total number of bonds in s */
256 sortable *s, /* bidirectional list of bonds */
257 t_nextnb *nnb) /* the tmp storage for excl */
258 /* Assume excl is initalised and s[] contains all bonds bidirectional */
263 for (i = 0; (i < nnb->nr); i++)
265 add_nnb(nnb, 0, i, i);
267 print_nnb(nnb, "After exclude self");
269 /* exclude all the bonded atoms */
272 for (i = 0; (i < nrbonds); i++)
274 add_nnb(nnb, 1, s[i].ai, s[i].aj);
277 print_nnb(nnb, "After exclude bonds");
279 /* for the nr of exclusions per atom */
280 for (n = 1; (n < nnb->nrex); n++)
282 /* now for all atoms */
283 for (i = 0; (i < nnb->nr); i++)
285 /* for all directly bonded atoms of atom i */
286 for (j = 0; (j < nnb->nrexcl[i][1]); j++)
289 /* store the 1st neighbour in nb */
290 nb = nnb->a[i][1][j];
292 /* store all atoms in nb's n-th list into i's n+1-th list */
293 for (k = 0; (k < nnb->nrexcl[nb][n]); k++)
295 if (i != nnb->a[nb][n][k])
297 add_nnb(nnb, n+1, i, nnb->a[nb][n][k]);
303 print_nnb(nnb, "After exclude rest");
307 static void add_b(t_params *bonds, int *nrf, sortable *s)
312 for (i = 0; (i < bonds->nr); i++)
314 ai = bonds->param[i].AI;
315 aj = bonds->param[i].AJ;
316 if ((ai < 0) || (aj < 0))
318 gmx_fatal(FARGS, "Impossible atom numbers in bond %d: ai=%d, aj=%d",
321 /* Add every bond twice */
329 void gen_nnb(t_nextnb *nnb, t_params plist[])
335 for (i = 0; (i < F_NRE); i++)
339 /* we need every bond twice (bidirectional) */
340 nrbonds += 2*plist[i].nr;
347 for (i = 0; (i < F_NRE); i++)
351 add_b(&plist[i], &nrf, s);
355 /* now sort the bonds */
356 prints("gen_excl before qsort", nrbonds, s);
359 qsort((void *) s, nrbonds, (size_t)sizeof(sortable), bond_sort);
360 prints("gen_excl after qsort", nrbonds, s);
363 do_gen(nrbonds, s, nnb);
368 sort_and_purge_nnb(t_nextnb *nnb)
370 int i, j, k, m, n, cnt, found, prev, idx;
372 for (i = 0; (i < nnb->nr); i++)
374 for (n = 0; (n <= nnb->nrex); n++)
376 /* Sort atoms in this list */
377 qsort(nnb->a[i][n], nnb->nrexcl[i][n], sizeof(int), compare_int);
381 for (j = 0; j < nnb->nrexcl[i][n]; j++)
383 idx = nnb->a[i][n][j];
386 for (m = 0; m < n && !found; m++)
388 for (k = 0; k < nnb->nrexcl[i][m] && !found; k++)
390 found = (idx == nnb->a[i][m][k]);
394 if (!found && nnb->a[i][n][j] != prev)
396 nnb->a[i][n][cnt] = nnb->a[i][n][j];
397 prev = nnb->a[i][n][cnt];
401 nnb->nrexcl[i][n] = cnt;
407 void generate_excl (int nrexcl, int nratoms, t_params plist[], t_nextnb *nnb, t_blocka *excl)
413 gmx_fatal(FARGS, "Can't have %d exclusions...", nrexcl);
415 init_nnb(nnb, nratoms, nrexcl);
418 sort_and_purge_nnb(nnb);
419 nnb2excl (nnb, excl);