Fix undefined behavior flagged by UBSAN
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / gpp_nextnb.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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,2020, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 /* This file is completely threadsafe - keep it that way! */
38 #include "gmxpre.h"
39
40 #include "gpp_nextnb.h"
41
42 #include <cstdlib>
43
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"
50
51 /* #define DEBUG_NNB */
52
53 typedef struct
54 {
55     int ai, aj;
56 } sortable;
57
58 static int bond_sort(const void* a, const void* b)
59 {
60     const sortable *sa, *sb;
61
62     sa = reinterpret_cast<const sortable*>(a);
63     sb = reinterpret_cast<const sortable*>(b);
64
65     if (sa->ai == sb->ai)
66     {
67         return (sa->aj - sb->aj);
68     }
69     else
70     {
71         return (sa->ai - sb->ai);
72     }
73 }
74
75 static int compare_int(const void* a, const void* b)
76 {
77     return (*reinterpret_cast<const int*>(a) - *reinterpret_cast<const int*>(b));
78 }
79
80
81 #ifdef DEBUG
82 #    define prints(str, n, s) __prints(str, n, s)
83 static void __prints(char* str, int n, sortable* s)
84 {
85     int i;
86
87     if (debug)
88     {
89         fprintf(debug, "%s\n", str);
90         fprintf(debug, "Sortables \n");
91         for (i = 0; (i < n); i++)
92         {
93             fprintf(debug, "%d\t%d\n", s[i].ai, s[i].aj);
94         }
95
96         fflush(debug);
97     }
98 }
99 #else
100 #    define prints(str, n, s)
101 #endif
102
103 void init_nnb(t_nextnb* nnb, int nr, int nrex)
104 {
105     int i;
106
107     /* initiate nnb */
108     nnb->nr   = nr;
109     nnb->nrex = nrex;
110
111     snew(nnb->a, nr);
112     snew(nnb->nrexcl, nr);
113     for (i = 0; (i < nr); i++)
114     {
115         snew(nnb->a[i], nrex + 1);
116         snew(nnb->nrexcl[i], nrex + 1);
117     }
118 }
119
120 static void add_nnb(t_nextnb* nnb, int nre, int i, int j)
121 {
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]++;
125 }
126
127 void done_nnb(t_nextnb* nnb)
128 {
129     int i, nre;
130
131     for (i = 0; (i < nnb->nr); i++)
132     {
133         for (nre = 0; (nre <= nnb->nrex); nre++)
134         {
135             if (nnb->nrexcl[i][nre] > 0)
136             {
137                 sfree(nnb->a[i][nre]);
138             }
139         }
140         sfree(nnb->nrexcl[i]);
141         sfree(nnb->a[i]);
142     }
143     sfree(nnb->a);
144     sfree(nnb->nrexcl);
145     nnb->nr   = 0;
146     nnb->nrex = 0;
147 }
148
149 #ifdef DEBUG_NNB
150 void __print_nnb(t_nextnb* nnb, char* s)
151 {
152     int i, j, k;
153
154     if (debug)
155     {
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++)
160         {
161             for (j = 0; (j <= nnb->nrex); j++)
162             {
163                 fprintf(debug, "nrexcl[%d][%d]: %d, excl: ", i, j, nnb->nrexcl[i][j]);
164                 for (k = 0; (k < nnb->nrexcl[i][j]); k++)
165                 {
166                     fprintf(debug, "%d, ", nnb->a[i][j][k]);
167                 }
168                 fprintf(debug, "\n");
169             }
170         }
171     }
172 }
173 #endif
174
175 static void nnb2excl(t_nextnb* nnb, gmx::ListOfLists<int>* excls)
176 {
177     int       i, j, j_index;
178     int       nre, nrx, nrs, nr_of_sortables;
179     sortable* s;
180
181     excls->clear();
182
183     for (i = 0; (i < nnb->nr); i++)
184     {
185         /* calculate the total number of exclusions for atom i */
186         nr_of_sortables = 0;
187         for (nre = 0; (nre <= nnb->nrex); nre++)
188         {
189             nr_of_sortables += nnb->nrexcl[i][nre];
190         }
191
192         /* make space for sortable array */
193         snew(s, nr_of_sortables);
194
195         /* fill the sortable array and sort it */
196         nrs = 0;
197         for (nre = 0; (nre <= nnb->nrex); nre++)
198         {
199             for (nrx = 0; (nrx < nnb->nrexcl[i][nre]); nrx++)
200             {
201                 s[nrs].ai = i;
202                 s[nrs].aj = nnb->a[i][nre][nrx];
203                 nrs++;
204             }
205         }
206         if (nrs != nr_of_sortables)
207         {
208             gmx_incons("Generating exclusions");
209         }
210         prints("nnb2excl before qsort", nr_of_sortables, s);
211         if (nr_of_sortables > 1)
212         {
213             qsort(s, nr_of_sortables, static_cast<size_t>(sizeof(s[0])), bond_sort);
214             prints("nnb2excl after qsort", nr_of_sortables, s);
215         }
216
217         /* remove duplicate entries from the list */
218         j_index = 0;
219         if (nr_of_sortables > 0)
220         {
221             for (j = 1; (j < nr_of_sortables); j++)
222             {
223                 if ((s[j].ai != s[j - 1].ai) || (s[j].aj != s[j - 1].aj))
224                 {
225                     s[j_index++] = s[j - 1];
226                 }
227             }
228             s[j_index++] = s[j - 1];
229         }
230         nr_of_sortables = j_index;
231         prints("after rm-double", j_index, s);
232
233         /* put the sorted exclusions in the target list */
234         excls->pushBackListOfSize(nr_of_sortables);
235         gmx::ArrayRef<int> exclusionsForAtom = excls->back();
236         for (nrs = 0; (nrs < nr_of_sortables); nrs++)
237         {
238             exclusionsForAtom[nrs] = s[nrs].aj;
239         }
240
241         /* cleanup temporary space */
242         sfree(s);
243     }
244 }
245
246 /*! \brief Return true of neighbor is already present in some exclusion level
247  *
248  * To avoid exploding complexity when processing exclusions for highly
249  * connected molecules with lots of exclusions, this routine is used to
250  * check whether a particular neighbor has already been excluded at any lower
251  * bond distance, in which case we should not add it to avoid creating loops.
252  *
253  * \param nnb            Valid initialized next-neighbor structure
254  * \param atom           The host atom whose neighbors we are searching
255  * \param highest_order  The highest-rank neighbor list to search.
256  * \param query          Atom index to look for
257  *
258  * \return True if query is present as an exclusion of up to highest_order
259  *         (inclusive) from atom. For instance, if highest_order is 2,
260  *         the routine will return true if the query atom is already listed as
261  *         first or second neighbor (exclusion) in nnb.
262  */
263 static bool atom_is_present_in_nnb(const t_nextnb* nnb, int atom, int highest_order, int query)
264 {
265     GMX_RELEASE_ASSERT(highest_order < nnb->nrex, "Inconsistent nnb seach parameters");
266
267     for (int order = 0; order <= highest_order; order++)
268     {
269         for (int m = 0; m < nnb->nrexcl[atom][order]; m++)
270         {
271             if (nnb->a[atom][order][m] == query)
272             {
273                 return true;
274             }
275         }
276     }
277     return false;
278 }
279
280 static void do_gen(int       nrbonds, /* total number of bonds in s     */
281                    sortable* s,       /* bidirectional list of bonds    */
282                    t_nextnb* nnb)     /* the tmp storage for excl     */
283 /* Assume excl is initalised and s[] contains all bonds bidirectional */
284 {
285     int i, j, k, n, nb;
286
287     /* exclude self */
288     for (i = 0; (i < nnb->nr); i++)
289     {
290         add_nnb(nnb, 0, i, i);
291     }
292     print_nnb(nnb, "After exclude self");
293
294     /* exclude all the bonded atoms */
295     if (nnb->nrex > 0)
296     {
297         for (i = 0; (i < nrbonds); i++)
298         {
299             add_nnb(nnb, 1, s[i].ai, s[i].aj);
300         }
301     }
302     print_nnb(nnb, "After exclude bonds");
303
304     /* for the nr of exclusions per atom */
305     for (n = 1; (n < nnb->nrex); n++)
306     {
307         /* now for all atoms */
308         for (i = 0; (i < nnb->nr); i++)
309         {
310             /* for all directly bonded atoms of atom i */
311             for (j = 0; (j < nnb->nrexcl[i][1]); j++)
312             {
313
314                 /* store the 1st neighbour in nb */
315                 nb = nnb->a[i][1][j];
316
317                 /* store all atoms in nb's n-th list into i's n+1-th list */
318                 for (k = 0; (k < nnb->nrexcl[nb][n]); k++)
319                 {
320                     // Only add if it is not already present as a closer neighbor
321                     // to avoid exploding complexity for highly connected molecules
322                     // with high exclusion order
323                     if (!atom_is_present_in_nnb(nnb, i, n, nnb->a[nb][n][k]))
324                     {
325                         add_nnb(nnb, n + 1, i, nnb->a[nb][n][k]);
326                     }
327                 }
328             }
329         }
330     }
331     print_nnb(nnb, "After exclude rest");
332 }
333
334 static void add_b(InteractionsOfType* bonds, int* nrf, sortable* s)
335 {
336     int i = 0;
337     for (const auto& bond : bonds->interactionTypes)
338     {
339         int ai = bond.ai();
340         int aj = bond.aj();
341         if ((ai < 0) || (aj < 0))
342         {
343             gmx_fatal(FARGS, "Impossible atom numbers in bond %d: ai=%d, aj=%d", i, ai, aj);
344         }
345         /* Add every bond twice */
346         s[(*nrf)].ai   = ai;
347         s[(*nrf)++].aj = aj;
348         s[(*nrf)].aj   = ai;
349         s[(*nrf)++].ai = aj;
350         i++;
351     }
352 }
353
354 void gen_nnb(t_nextnb* nnb, gmx::ArrayRef<InteractionsOfType> plist)
355 {
356     sortable* s;
357     int       nrbonds, nrf;
358
359     nrbonds = 0;
360     for (int i = 0; (i < F_NRE); i++)
361     {
362         if (IS_CHEMBOND(i))
363         {
364             /* we need every bond twice (bidirectional) */
365             nrbonds += 2 * plist[i].size();
366         }
367     }
368
369     snew(s, nrbonds);
370
371     nrf = 0;
372     for (int i = 0; (i < F_NRE); i++)
373     {
374         if (IS_CHEMBOND(i))
375         {
376             add_b(&plist[i], &nrf, s);
377         }
378     }
379
380     /* now sort the bonds */
381     prints("gen_excl before qsort", nrbonds, s);
382     if (nrbonds > 1)
383     {
384         qsort(s, nrbonds, static_cast<size_t>(sizeof(sortable)), bond_sort);
385         prints("gen_excl after qsort", nrbonds, s);
386     }
387
388     do_gen(nrbonds, s, nnb);
389     sfree(s);
390 }
391
392 static void sort_and_purge_nnb(t_nextnb* nnb)
393 {
394     int  i, j, k, m, n, cnt, prev, idx;
395     bool found;
396
397     for (i = 0; (i < nnb->nr); i++)
398     {
399         for (n = 0; (n <= nnb->nrex); n++)
400         {
401             /* Sort atoms in this list */
402             if (nnb->nrexcl[i][n] > 0)
403             {
404                 qsort(nnb->a[i][n], nnb->nrexcl[i][n], sizeof(int), compare_int);
405             }
406             cnt  = 0;
407             prev = -1;
408             for (j = 0; j < nnb->nrexcl[i][n]; j++)
409             {
410                 idx = nnb->a[i][n][j];
411
412                 found = false;
413                 for (m = 0; m < n && !found; m++)
414                 {
415                     for (k = 0; k < nnb->nrexcl[i][m] && !found; k++)
416                     {
417                         found = idx == nnb->a[i][m][k];
418                     }
419                 }
420
421                 if (!found && nnb->a[i][n][j] != prev)
422                 {
423                     nnb->a[i][n][cnt] = nnb->a[i][n][j];
424                     prev              = nnb->a[i][n][cnt];
425                     cnt++;
426                 }
427             }
428             nnb->nrexcl[i][n] = cnt;
429         }
430     }
431 }
432
433
434 void generate_excl(int nrexcl, int nratoms, gmx::ArrayRef<InteractionsOfType> plist, gmx::ListOfLists<int>* excls)
435 {
436     t_nextnb nnb;
437     if (nrexcl < 0)
438     {
439         gmx_fatal(FARGS, "Can't have %d exclusions...", nrexcl);
440     }
441     init_nnb(&nnb, nratoms, nrexcl);
442     gen_nnb(&nnb, plist);
443     sort_and_purge_nnb(&nnb);
444     nnb2excl(&nnb, excls);
445     done_nnb(&nnb);
446 }