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 static int 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);
75 static int compare_int(const void* a, const void* b)
77 return (*reinterpret_cast<const int*>(a) - *reinterpret_cast<const 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];
192 /* make space for sortable array */
193 snew(s, nr_of_sortables);
195 /* fill the sortable array and sort it */
197 for (nre = 0; (nre <= nnb->nrex); nre++)
199 for (nrx = 0; (nrx < nnb->nrexcl[i][nre]); nrx++)
202 s[nrs].aj = nnb->a[i][nre][nrx];
206 if (nrs != nr_of_sortables)
208 gmx_incons("Generating exclusions");
210 prints("nnb2excl before qsort", nr_of_sortables, s);
211 if (nr_of_sortables > 1)
213 qsort(s, nr_of_sortables, static_cast<size_t>(sizeof(s[0])), bond_sort);
214 prints("nnb2excl after qsort", nr_of_sortables, s);
217 /* remove duplicate entries from the list */
219 if (nr_of_sortables > 0)
221 for (j = 1; (j < nr_of_sortables); j++)
223 if ((s[j].ai != s[j - 1].ai) || (s[j].aj != s[j - 1].aj))
225 s[j_index++] = s[j - 1];
228 s[j_index++] = s[j - 1];
230 nr_of_sortables = j_index;
231 prints("after rm-double", j_index, s);
233 /* make space for arrays */
234 srenew(excl->a, excl->nra + nr_of_sortables);
236 /* put the sorted exclusions in the target list */
237 for (nrs = 0; (nrs < nr_of_sortables); nrs++)
239 excl->a[excl->nra + nrs] = s[nrs].aj;
241 excl->nra += nr_of_sortables;
242 excl->index[i + 1] = excl->nra;
244 /* cleanup temporary space */
249 /*! \brief Return true of neighbor is already present in some exclusion level
251 * To avoid exploding complexity when processing exclusions for highly
252 * connected molecules with lots of exclusions, this routine is used to
253 * check whether a particular neighbor has already been excluded at any lower
254 * bond distance, in which case we should not add it to avoid creating loops.
256 * \param nnb Valid initialized next-neighbor structure
257 * \param atom The host atom whose neighbors we are searching
258 * \param highest_order The highest-rank neighbor list to search.
259 * \param query Atom index to look for
261 * \return True if query is present as an exclusion of up to highest_order
262 * (inclusive) from atom. For instance, if highest_order is 2,
263 * the routine will return true if the query atom is already listed as
264 * first or second neighbor (exclusion) in nnb.
266 static bool atom_is_present_in_nnb(const t_nextnb* nnb, int atom, int highest_order, int query)
268 GMX_RELEASE_ASSERT(highest_order < nnb->nrex, "Inconsistent nnb seach parameters");
270 for (int order = 0; order <= highest_order; order++)
272 for (int m = 0; m < nnb->nrexcl[atom][order]; m++)
274 if (nnb->a[atom][order][m] == query)
283 static void do_gen(int nrbonds, /* total number of bonds in s */
284 sortable* s, /* bidirectional list of bonds */
285 t_nextnb* nnb) /* the tmp storage for excl */
286 /* Assume excl is initalised and s[] contains all bonds bidirectional */
291 for (i = 0; (i < nnb->nr); i++)
293 add_nnb(nnb, 0, i, i);
295 print_nnb(nnb, "After exclude self");
297 /* exclude all the bonded atoms */
300 for (i = 0; (i < nrbonds); i++)
302 add_nnb(nnb, 1, s[i].ai, s[i].aj);
305 print_nnb(nnb, "After exclude bonds");
307 /* for the nr of exclusions per atom */
308 for (n = 1; (n < nnb->nrex); n++)
310 /* now for all atoms */
311 for (i = 0; (i < nnb->nr); i++)
313 /* for all directly bonded atoms of atom i */
314 for (j = 0; (j < nnb->nrexcl[i][1]); j++)
317 /* store the 1st neighbour in nb */
318 nb = nnb->a[i][1][j];
320 /* store all atoms in nb's n-th list into i's n+1-th list */
321 for (k = 0; (k < nnb->nrexcl[nb][n]); k++)
323 // Only add if it is not already present as a closer neighbor
324 // to avoid exploding complexity for highly connected molecules
325 // with high exclusion order
326 if (!atom_is_present_in_nnb(nnb, i, n, nnb->a[nb][n][k]))
328 add_nnb(nnb, n + 1, i, nnb->a[nb][n][k]);
334 print_nnb(nnb, "After exclude rest");
337 static void add_b(InteractionsOfType* bonds, int* nrf, sortable* s)
340 for (const auto& bond : bonds->interactionTypes)
344 if ((ai < 0) || (aj < 0))
346 gmx_fatal(FARGS, "Impossible atom numbers in bond %d: ai=%d, aj=%d", i, ai, aj);
348 /* Add every bond twice */
357 void gen_nnb(t_nextnb* nnb, gmx::ArrayRef<InteractionsOfType> plist)
363 for (int i = 0; (i < F_NRE); i++)
367 /* we need every bond twice (bidirectional) */
368 nrbonds += 2 * plist[i].size();
375 for (int i = 0; (i < F_NRE); i++)
379 add_b(&plist[i], &nrf, s);
383 /* now sort the bonds */
384 prints("gen_excl before qsort", nrbonds, s);
387 qsort(s, nrbonds, static_cast<size_t>(sizeof(sortable)), bond_sort);
388 prints("gen_excl after qsort", nrbonds, s);
391 do_gen(nrbonds, s, nnb);
395 static void sort_and_purge_nnb(t_nextnb* nnb)
397 int i, j, k, m, n, cnt, prev, idx;
400 for (i = 0; (i < nnb->nr); i++)
402 for (n = 0; (n <= nnb->nrex); n++)
404 /* Sort atoms in this list */
405 qsort(nnb->a[i][n], nnb->nrexcl[i][n], sizeof(int), compare_int);
409 for (j = 0; j < nnb->nrexcl[i][n]; j++)
411 idx = nnb->a[i][n][j];
414 for (m = 0; m < n && !found; m++)
416 for (k = 0; k < nnb->nrexcl[i][m] && !found; k++)
418 found = idx == nnb->a[i][m][k];
422 if (!found && nnb->a[i][n][j] != prev)
424 nnb->a[i][n][cnt] = nnb->a[i][n][j];
425 prev = nnb->a[i][n][cnt];
429 nnb->nrexcl[i][n] = cnt;
435 void generate_excl(int nrexcl, int nratoms, gmx::ArrayRef<InteractionsOfType> plist, t_blocka* excl)
440 gmx_fatal(FARGS, "Can't have %d exclusions...", nrexcl);
442 init_nnb(&nnb, nratoms, nrexcl);
443 gen_nnb(&nnb, plist);
445 sort_and_purge_nnb(&nnb);
446 nnb2excl(&nnb, excl);