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,2015,2018,2019, 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! */
40 #include "gpp_nextnb.h"
44 #include "gromacs/gmxpreprocess/grompp_impl.h"
45 #include "gromacs/gmxpreprocess/toputil.h"
46 #include "gromacs/topology/ifunc.h"
47 #include "gromacs/utility/fatalerror.h"
48 #include "gromacs/utility/gmxassert.h"
49 #include "gromacs/utility/smalloc.h"
51 /* #define DEBUG_NNB */
58 bond_sort (const void *a, const void *b)
60 const sortable *sa, *sb;
62 sa = reinterpret_cast<const sortable*>(a);
63 sb = reinterpret_cast<const sortable*>(b);
67 return (sa->aj-sb->aj);
71 return (sa->ai-sb->ai);
76 compare_int (const void * a, const void * b)
78 return ( *reinterpret_cast<const int*>(a) - *reinterpret_cast<const int*>(b) );
83 #define prints(str, n, s) __prints(str, n, s)
84 static void __prints(char *str, int n, sortable *s)
90 fprintf(debug, "%s\n", str);
91 fprintf(debug, "Sortables \n");
92 for (i = 0; (i < n); i++)
94 fprintf(debug, "%d\t%d\n", s[i].ai, s[i].aj);
101 #define prints(str, n, s)
104 void init_nnb(t_nextnb *nnb, int nr, int nrex)
113 snew(nnb->nrexcl, nr);
114 for (i = 0; (i < nr); i++)
116 snew(nnb->a[i], nrex+1);
117 snew(nnb->nrexcl[i], nrex+1);
121 static void add_nnb (t_nextnb *nnb, int nre, int i, int j)
123 srenew(nnb->a[i][nre], nnb->nrexcl[i][nre]+1);
124 nnb->a[i][nre][nnb->nrexcl[i][nre]] = j;
125 nnb->nrexcl[i][nre]++;
128 void done_nnb (t_nextnb *nnb)
132 for (i = 0; (i < nnb->nr); i++)
134 for (nre = 0; (nre <= nnb->nrex); nre++)
136 if (nnb->nrexcl[i][nre] > 0)
138 sfree (nnb->a[i][nre]);
141 sfree (nnb->nrexcl[i]);
151 void __print_nnb(t_nextnb *nnb, char *s)
157 fprintf(debug, "%s\n", s);
158 fprintf(debug, "nnb->nr: %d\n", nnb->nr);
159 fprintf(debug, "nnb->nrex: %d\n", nnb->nrex);
160 for (i = 0; (i < nnb->nr); i++)
162 for (j = 0; (j <= nnb->nrex); j++)
164 fprintf(debug, "nrexcl[%d][%d]: %d, excl: ", i, j, nnb->nrexcl[i][j]);
165 for (k = 0; (k < nnb->nrexcl[i][j]); k++)
167 fprintf(debug, "%d, ", nnb->a[i][j][k]);
169 fprintf(debug, "\n");
176 static void nnb2excl(t_nextnb *nnb, t_blocka *excl)
179 int nre, nrx, nrs, nr_of_sortables;
182 srenew(excl->index, nnb->nr+1);
184 for (i = 0; (i < nnb->nr); i++)
186 /* calculate the total number of exclusions for atom i */
188 for (nre = 0; (nre <= nnb->nrex); nre++)
190 nr_of_sortables += nnb->nrexcl[i][nre];
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 (s, nr_of_sortables, static_cast<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 /*! \brief Return true of neighbor is already present in some exclusion level
252 * To avoid exploding complexity when processing exclusions for highly
253 * connected molecules with lots of exclusions, this routine is used to
254 * check whether a particular neighbor has already been excluded at any lower
255 * bond distance, in which case we should not add it to avoid creating loops.
257 * \param nnb Valid initialized next-neighbor structure
258 * \param atom The host atom whose neighbors we are searching
259 * \param highest_order The highest-rank neighbor list to search.
260 * \param query Atom index to look for
262 * \return True if query is present as an exclusion of up to highest_order
263 * (inclusive) from atom. For instance, if highest_order is 2,
264 * the routine will return true if the query atom is already listed as
265 * first or second neighbor (exclusion) in nnb.
268 atom_is_present_in_nnb(const t_nextnb * nnb,
273 GMX_RELEASE_ASSERT(highest_order < nnb->nrex, "Inconsistent nnb seach parameters");
275 for (int order = 0; order <= highest_order; order++)
277 for (int m = 0; m < nnb->nrexcl[atom][order]; m++)
279 if (nnb->a[atom][order][m] == query)
288 static void do_gen(int nrbonds, /* total number of bonds in s */
289 sortable *s, /* bidirectional list of bonds */
290 t_nextnb *nnb) /* the tmp storage for excl */
291 /* Assume excl is initalised and s[] contains all bonds bidirectional */
296 for (i = 0; (i < nnb->nr); i++)
298 add_nnb(nnb, 0, i, i);
300 print_nnb(nnb, "After exclude self");
302 /* exclude all the bonded atoms */
305 for (i = 0; (i < nrbonds); i++)
307 add_nnb(nnb, 1, s[i].ai, s[i].aj);
310 print_nnb(nnb, "After exclude bonds");
312 /* for the nr of exclusions per atom */
313 for (n = 1; (n < nnb->nrex); n++)
315 /* now for all atoms */
316 for (i = 0; (i < nnb->nr); i++)
318 /* for all directly bonded atoms of atom i */
319 for (j = 0; (j < nnb->nrexcl[i][1]); j++)
322 /* store the 1st neighbour in nb */
323 nb = nnb->a[i][1][j];
325 /* store all atoms in nb's n-th list into i's n+1-th list */
326 for (k = 0; (k < nnb->nrexcl[nb][n]); k++)
328 // Only add if it is not already present as a closer neighbor
329 // to avoid exploding complexity for highly connected molecules
330 // with high exclusion order
331 if (!atom_is_present_in_nnb(nnb, i, n, nnb->a[nb][n][k]))
333 add_nnb(nnb, n+1, i, nnb->a[nb][n][k]);
339 print_nnb(nnb, "After exclude rest");
343 static void add_b(t_params *bonds, int *nrf, sortable *s)
348 for (i = 0; (i < bonds->nr); i++)
350 ai = bonds->param[i].ai();
351 aj = bonds->param[i].aj();
352 if ((ai < 0) || (aj < 0))
354 gmx_fatal(FARGS, "Impossible atom numbers in bond %d: ai=%d, aj=%d",
357 /* Add every bond twice */
365 void gen_nnb(t_nextnb *nnb, t_params plist[])
371 for (i = 0; (i < F_NRE); i++)
375 /* we need every bond twice (bidirectional) */
376 nrbonds += 2*plist[i].nr;
383 for (i = 0; (i < F_NRE); i++)
387 add_b(&plist[i], &nrf, s);
391 /* now sort the bonds */
392 prints("gen_excl before qsort", nrbonds, s);
395 qsort(s, nrbonds, static_cast<size_t>(sizeof(sortable)), bond_sort);
396 prints("gen_excl after qsort", nrbonds, s);
399 do_gen(nrbonds, s, nnb);
404 sort_and_purge_nnb(t_nextnb *nnb)
406 int i, j, k, m, n, cnt, prev, idx;
409 for (i = 0; (i < nnb->nr); i++)
411 for (n = 0; (n <= nnb->nrex); n++)
413 /* Sort atoms in this list */
414 qsort(nnb->a[i][n], nnb->nrexcl[i][n], sizeof(int), compare_int);
418 for (j = 0; j < nnb->nrexcl[i][n]; j++)
420 idx = nnb->a[i][n][j];
423 for (m = 0; m < n && !found; m++)
425 for (k = 0; k < nnb->nrexcl[i][m] && !found; k++)
427 found = idx == nnb->a[i][m][k];
431 if (!found && nnb->a[i][n][j] != prev)
433 nnb->a[i][n][cnt] = nnb->a[i][n][j];
434 prev = nnb->a[i][n][cnt];
438 nnb->nrexcl[i][n] = cnt;
444 void generate_excl (int nrexcl, int nratoms, t_params plist[], t_nextnb *nnb, t_blocka *excl)
448 gmx_fatal(FARGS, "Can't have %d exclusions...", nrexcl);
450 init_nnb(nnb, nratoms, nrexcl);
453 sort_and_purge_nnb(nnb);
454 nnb2excl (nnb, excl);