Fix part of old-style casting
[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, 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 <stdlib.h>
43
44 #include "gromacs/gmxpreprocess/toputil.h"
45 #include "gromacs/topology/ifunc.h"
46 #include "gromacs/utility/fatalerror.h"
47 #include "gromacs/utility/gmxassert.h"
48 #include "gromacs/utility/smalloc.h"
49
50 /* #define DEBUG_NNB */
51
52 typedef struct {
53     int ai, aj;
54 } sortable;
55
56 static int
57 bond_sort (const void *a, const void *b)
58 {
59     sortable *sa, *sb;
60
61     sa = (sortable *) a;
62     sb = (sortable *) b;
63
64     if (sa->ai == sb->ai)
65     {
66         return (sa->aj-sb->aj);
67     }
68     else
69     {
70         return (sa->ai-sb->ai);
71     }
72 }
73
74 static int
75 compare_int (const void * a, const void * b)
76 {
77     return ( *(int*)a - *(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, t_blocka *excl)
176 {
177     int       i, j, j_index;
178     int       nre, nrx, nrs, nr_of_sortables;
179     sortable *s;
180
181     srenew(excl->index, nnb->nr+1);
182     excl->index[0] = 0;
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 ((void *)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         /* make space for arrays */
234         srenew(excl->a, excl->nra+nr_of_sortables);
235
236         /* put the sorted exclusions in the target list */
237         for (nrs = 0; (nrs < nr_of_sortables); nrs++)
238         {
239             excl->a[excl->nra+nrs] = s[nrs].aj;
240         }
241         excl->nra       += nr_of_sortables;
242         excl->index[i+1] = excl->nra;
243
244         /* cleanup temporary space */
245         sfree (s);
246     }
247 }
248
249 /*! \brief Return true of neighbor is already present in some exclusion level
250  *
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.
255  *
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
260  *
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.
265  */
266 static bool
267 atom_is_present_in_nnb(const t_nextnb *   nnb,
268                        int                atom,
269                        int                highest_order,
270                        int                query)
271 {
272     GMX_RELEASE_ASSERT(highest_order < nnb->nrex, "Inconsistent nnb seach parameters");
273
274     for (int order = 0; order <= highest_order; order++)
275     {
276         for (int m = 0; m < nnb->nrexcl[atom][order]; m++)
277         {
278             if (nnb->a[atom][order][m] == query)
279             {
280                 return true;
281             }
282         }
283     }
284     return false;
285 }
286
287 static void do_gen(int       nrbonds, /* total number of bonds in s     */
288                    sortable *s,       /* bidirectional list of bonds    */
289                    t_nextnb *nnb)     /* the tmp storage for excl     */
290 /* Assume excl is initalised and s[] contains all bonds bidirectional */
291 {
292     int i, j, k, n, nb;
293
294     /* exclude self */
295     for (i = 0; (i < nnb->nr); i++)
296     {
297         add_nnb(nnb, 0, i, i);
298     }
299     print_nnb(nnb, "After exclude self");
300
301     /* exclude all the bonded atoms */
302     if (nnb->nrex > 0)
303     {
304         for (i = 0; (i < nrbonds); i++)
305         {
306             add_nnb(nnb, 1, s[i].ai, s[i].aj);
307         }
308     }
309     print_nnb(nnb, "After exclude bonds");
310
311     /* for the nr of exclusions per atom */
312     for (n = 1; (n < nnb->nrex); n++)
313     {
314         /* now for all atoms */
315         for (i = 0; (i < nnb->nr); i++)
316         {
317             /* for all directly bonded atoms of atom i */
318             for (j = 0; (j < nnb->nrexcl[i][1]); j++)
319             {
320
321                 /* store the 1st neighbour in nb */
322                 nb = nnb->a[i][1][j];
323
324                 /* store all atoms in nb's n-th list into i's n+1-th list */
325                 for (k = 0; (k < nnb->nrexcl[nb][n]); k++)
326                 {
327                     // Only add if it is not already present as a closer neighbor
328                     // to avoid exploding complexity for highly connected molecules
329                     // with high exclusion order
330                     if (!atom_is_present_in_nnb(nnb, i, n, nnb->a[nb][n][k]))
331                     {
332                         add_nnb(nnb, n+1, i, nnb->a[nb][n][k]);
333                     }
334                 }
335             }
336         }
337     }
338     print_nnb(nnb, "After exclude rest");
339
340 }
341
342 static void add_b(t_params *bonds, int *nrf, sortable *s)
343 {
344     int i;
345     int ai, aj;
346
347     for (i = 0; (i < bonds->nr); i++)
348     {
349         ai = bonds->param[i].ai();
350         aj = bonds->param[i].aj();
351         if ((ai < 0) || (aj < 0))
352         {
353             gmx_fatal(FARGS, "Impossible atom numbers in bond %d: ai=%d, aj=%d",
354                       i, ai, aj);
355         }
356         /* Add every bond twice */
357         s[(*nrf)].ai   = ai;
358         s[(*nrf)++].aj = aj;
359         s[(*nrf)].aj   = ai;
360         s[(*nrf)++].ai = aj;
361     }
362 }
363
364 void gen_nnb(t_nextnb *nnb, t_params plist[])
365 {
366     sortable *s;
367     int       i, nrbonds, nrf;
368
369     nrbonds = 0;
370     for (i = 0; (i < F_NRE); i++)
371     {
372         if (IS_CHEMBOND(i))
373         {
374             /* we need every bond twice (bidirectional) */
375             nrbonds += 2*plist[i].nr;
376         }
377     }
378
379     snew(s, nrbonds);
380
381     nrf = 0;
382     for (i = 0; (i < F_NRE); i++)
383     {
384         if (IS_CHEMBOND(i))
385         {
386             add_b(&plist[i], &nrf, s);
387         }
388     }
389
390     /* now sort the bonds */
391     prints("gen_excl before qsort", nrbonds, s);
392     if (nrbonds > 1)
393     {
394         qsort((void *) s, nrbonds, static_cast<size_t>(sizeof(sortable)), bond_sort);
395         prints("gen_excl after qsort", nrbonds, s);
396     }
397
398     do_gen(nrbonds, s, nnb);
399     sfree(s);
400 }
401
402 static void
403 sort_and_purge_nnb(t_nextnb *nnb)
404 {
405     int i, j, k, m, n, cnt, found, prev, idx;
406
407     for (i = 0; (i < nnb->nr); i++)
408     {
409         for (n = 0; (n <= nnb->nrex); n++)
410         {
411             /* Sort atoms in this list */
412             qsort(nnb->a[i][n], nnb->nrexcl[i][n], sizeof(int), compare_int);
413
414             cnt  = 0;
415             prev = -1;
416             for (j = 0; j < nnb->nrexcl[i][n]; j++)
417             {
418                 idx = nnb->a[i][n][j];
419
420                 found = 0;
421                 for (m = 0; m < n && !found; m++)
422                 {
423                     for (k = 0; k < nnb->nrexcl[i][m] && !found; k++)
424                     {
425                         found = (idx == nnb->a[i][m][k]);
426                     }
427                 }
428
429                 if (!found && nnb->a[i][n][j] != prev)
430                 {
431                     nnb->a[i][n][cnt] = nnb->a[i][n][j];
432                     prev              = nnb->a[i][n][cnt];
433                     cnt++;
434                 }
435             }
436             nnb->nrexcl[i][n] = cnt;
437         }
438     }
439 }
440
441
442 void generate_excl (int nrexcl, int nratoms, t_params plist[], t_nextnb *nnb, t_blocka *excl)
443 {
444     if (nrexcl < 0)
445     {
446         gmx_fatal(FARGS, "Can't have %d exclusions...", nrexcl);
447     }
448     init_nnb(nnb, nratoms, nrexcl);
449     gen_nnb(nnb, plist);
450     excl->nr = nratoms;
451     sort_and_purge_nnb(nnb);
452     nnb2excl (nnb, excl);
453 }