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"
49 #include "gromacs/utility/fatalerror.h"
50 #include "gromacs/utility/smalloc.h"
57 bond_sort (const void *a, const void *b)
66 return (sa->aj-sb->aj);
70 return (sa->ai-sb->ai);
75 compare_int (const void * a, const void * b)
77 return ( *(int*)a - *(int*)b );
82 #define prints(str, n, s) __prints(str, n, s)
83 static void __prints(char *str, int n, sortable *s)
89 fprintf(debug, "%s\n", str);
90 fprintf(debug, "Sortables \n");
91 for (i = 0; (i < n); i++)
93 fprintf(debug, "%d\t%d\n", s[i].ai, s[i].aj);
100 #define prints(str, n, s)
103 void init_nnb(t_nextnb *nnb, int nr, int nrex)
112 snew(nnb->nrexcl, nr);
113 for (i = 0; (i < nr); i++)
115 snew(nnb->a[i], nrex+1);
116 snew(nnb->nrexcl[i], nrex+1);
120 static void add_nnb (t_nextnb *nnb, int nre, int i, int j)
122 srenew(nnb->a[i][nre], nnb->nrexcl[i][nre]+1);
123 nnb->a[i][nre][nnb->nrexcl[i][nre]] = j;
124 nnb->nrexcl[i][nre]++;
127 void done_nnb (t_nextnb *nnb)
131 for (i = 0; (i < nnb->nr); i++)
133 for (nre = 0; (nre <= nnb->nrex); nre++)
135 if (nnb->nrexcl[i][nre] > 0)
137 sfree (nnb->a[i][nre]);
140 sfree (nnb->nrexcl[i]);
150 void __print_nnb(t_nextnb *nnb, char *s)
156 fprintf(debug, "%s\n", s);
157 fprintf(debug, "nnb->nr: %d\n", nnb->nr);
158 fprintf(debug, "nnb->nrex: %d\n", nnb->nrex);
159 for (i = 0; (i < nnb->nr); i++)
161 for (j = 0; (j <= nnb->nrex); j++)
163 fprintf(debug, "nrexcl[%d][%d]: %d, excl: ", i, j, nnb->nrexcl[i][j]);
164 for (k = 0; (k < nnb->nrexcl[i][j]); k++)
166 fprintf(debug, "%d, ", nnb->a[i][j][k]);
168 fprintf(debug, "\n");
175 static void nnb2excl(t_nextnb *nnb, t_blocka *excl)
178 int nre, nrx, nrs, nr_of_sortables;
181 srenew(excl->index, nnb->nr+1);
183 for (i = 0; (i < nnb->nr); i++)
185 /* calculate the total number of exclusions for atom i */
187 for (nre = 0; (nre <= nnb->nrex); nre++)
189 nr_of_sortables += nnb->nrexcl[i][nre];
194 fprintf(debug, "nr_of_sortables: %d\n", nr_of_sortables);
196 /* make space for sortable array */
197 snew(s, nr_of_sortables);
199 /* fill the sortable array and sort it */
201 for (nre = 0; (nre <= nnb->nrex); nre++)
203 for (nrx = 0; (nrx < nnb->nrexcl[i][nre]); nrx++)
206 s[nrs].aj = nnb->a[i][nre][nrx];
210 if (nrs != nr_of_sortables)
212 gmx_incons("Generating exclusions");
214 prints("nnb2excl before qsort", nr_of_sortables, s);
215 if (nr_of_sortables > 1)
217 qsort ((void *)s, nr_of_sortables, (size_t)sizeof(s[0]), bond_sort);
218 prints("nnb2excl after qsort", nr_of_sortables, s);
221 /* remove duplicate entries from the list */
223 if (nr_of_sortables > 0)
225 for (j = 1; (j < nr_of_sortables); j++)
227 if ((s[j].ai != s[j-1].ai) || (s[j].aj != s[j-1].aj))
229 s[j_index++] = s[j-1];
232 s[j_index++] = s[j-1];
234 nr_of_sortables = j_index;
235 prints("after rm-double", j_index, s);
237 /* make space for arrays */
238 srenew(excl->a, excl->nra+nr_of_sortables);
240 /* put the sorted exclusions in the target list */
241 for (nrs = 0; (nrs < nr_of_sortables); nrs++)
243 excl->a[excl->nra+nrs] = s[nrs].aj;
245 excl->nra += nr_of_sortables;
246 excl->index[i+1] = excl->nra;
248 /* cleanup temporary space */
253 print_blocka(debug, "Exclusions", "Atom", "Excluded", excl);
257 static void do_gen(int nrbonds, /* total number of bonds in s */
258 sortable *s, /* bidirectional list of bonds */
259 t_nextnb *nnb) /* the tmp storage for excl */
260 /* Assume excl is initalised and s[] contains all bonds bidirectional */
265 for (i = 0; (i < nnb->nr); i++)
267 add_nnb(nnb, 0, i, i);
269 print_nnb(nnb, "After exclude self");
271 /* exclude all the bonded atoms */
274 for (i = 0; (i < nrbonds); i++)
276 add_nnb(nnb, 1, s[i].ai, s[i].aj);
279 print_nnb(nnb, "After exclude bonds");
281 /* for the nr of exclusions per atom */
282 for (n = 1; (n < nnb->nrex); n++)
284 /* now for all atoms */
285 for (i = 0; (i < nnb->nr); i++)
287 /* for all directly bonded atoms of atom i */
288 for (j = 0; (j < nnb->nrexcl[i][1]); j++)
291 /* store the 1st neighbour in nb */
292 nb = nnb->a[i][1][j];
294 /* store all atoms in nb's n-th list into i's n+1-th list */
295 for (k = 0; (k < nnb->nrexcl[nb][n]); k++)
297 if (i != nnb->a[nb][n][k])
299 add_nnb(nnb, n+1, i, nnb->a[nb][n][k]);
305 print_nnb(nnb, "After exclude rest");
309 static void add_b(t_params *bonds, int *nrf, sortable *s)
314 for (i = 0; (i < bonds->nr); i++)
316 ai = bonds->param[i].AI;
317 aj = bonds->param[i].AJ;
318 if ((ai < 0) || (aj < 0))
320 gmx_fatal(FARGS, "Impossible atom numbers in bond %d: ai=%d, aj=%d",
323 /* Add every bond twice */
331 void gen_nnb(t_nextnb *nnb, t_params plist[])
337 for (i = 0; (i < F_NRE); i++)
341 /* we need every bond twice (bidirectional) */
342 nrbonds += 2*plist[i].nr;
349 for (i = 0; (i < F_NRE); i++)
353 add_b(&plist[i], &nrf, s);
357 /* now sort the bonds */
358 prints("gen_excl before qsort", nrbonds, s);
361 qsort((void *) s, nrbonds, (size_t)sizeof(sortable), bond_sort);
362 prints("gen_excl after qsort", nrbonds, s);
365 do_gen(nrbonds, s, nnb);
370 sort_and_purge_nnb(t_nextnb *nnb)
372 int i, j, k, m, n, cnt, found, prev, idx;
374 for (i = 0; (i < nnb->nr); i++)
376 for (n = 0; (n <= nnb->nrex); n++)
378 /* Sort atoms in this list */
379 qsort(nnb->a[i][n], nnb->nrexcl[i][n], sizeof(int), compare_int);
383 for (j = 0; j < nnb->nrexcl[i][n]; j++)
385 idx = nnb->a[i][n][j];
388 for (m = 0; m < n && !found; m++)
390 for (k = 0; k < nnb->nrexcl[i][m] && !found; k++)
392 found = (idx == nnb->a[i][m][k]);
396 if (!found && nnb->a[i][n][j] != prev)
398 nnb->a[i][n][cnt] = nnb->a[i][n][j];
399 prev = nnb->a[i][n][cnt];
403 nnb->nrexcl[i][n] = cnt;
409 void generate_excl (int nrexcl, int nratoms, t_params plist[], t_nextnb *nnb, t_blocka *excl)
415 gmx_fatal(FARGS, "Can't have %d exclusions...", nrexcl);
417 init_nnb(nnb, nratoms, nrexcl);
420 sort_and_purge_nnb(nnb);
421 nnb2excl (nnb, excl);