Manually sort some includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / gpp_nextnb.c
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, 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/utility/fatalerror.h"
46 #include "gromacs/utility/smalloc.h"
47
48 /* #define DEBUG_NNB */
49
50 typedef struct {
51     int ai, aj;
52 } sortable;
53
54 static int
55 bond_sort (const void *a, const void *b)
56 {
57     sortable *sa, *sb;
58
59     sa = (sortable *) a;
60     sb = (sortable *) b;
61
62     if (sa->ai == sb->ai)
63     {
64         return (sa->aj-sb->aj);
65     }
66     else
67     {
68         return (sa->ai-sb->ai);
69     }
70 }
71
72 static int
73 compare_int (const void * a, const void * b)
74 {
75     return ( *(int*)a - *(int*)b );
76 }
77
78
79 #ifdef DEBUG
80 #define prints(str, n, s) __prints(str, n, s)
81 static void __prints(char *str, int n, sortable *s)
82 {
83     int i;
84
85     if (debug)
86     {
87         fprintf(debug, "%s\n", str);
88         fprintf(debug, "Sortables \n");
89         for (i = 0; (i < n); i++)
90         {
91             fprintf(debug, "%d\t%d\n", s[i].ai, s[i].aj);
92         }
93
94         fflush(debug);
95     }
96 }
97 #else
98 #define prints(str, n, s)
99 #endif
100
101 void init_nnb(t_nextnb *nnb, int nr, int nrex)
102 {
103     int i;
104
105     /* initiate nnb */
106     nnb->nr   = nr;
107     nnb->nrex = nrex;
108
109     snew(nnb->a, nr);
110     snew(nnb->nrexcl, nr);
111     for (i = 0; (i < nr); i++)
112     {
113         snew(nnb->a[i], nrex+1);
114         snew(nnb->nrexcl[i], nrex+1);
115     }
116 }
117
118 static void add_nnb (t_nextnb *nnb, int nre, int i, int j)
119 {
120     srenew(nnb->a[i][nre], nnb->nrexcl[i][nre]+1);
121     nnb->a[i][nre][nnb->nrexcl[i][nre]] = j;
122     nnb->nrexcl[i][nre]++;
123 }
124
125 void done_nnb (t_nextnb *nnb)
126 {
127     int i, nre;
128
129     for (i = 0; (i < nnb->nr); i++)
130     {
131         for (nre = 0; (nre <= nnb->nrex); nre++)
132         {
133             if (nnb->nrexcl[i][nre] > 0)
134             {
135                 sfree (nnb->a[i][nre]);
136             }
137         }
138         sfree (nnb->nrexcl[i]);
139         sfree (nnb->a[i]);
140     }
141     sfree (nnb->a);
142     sfree (nnb->nrexcl);
143     nnb->nr   = 0;
144     nnb->nrex = 0;
145 }
146
147 #ifdef DEBUG_NNB
148 void __print_nnb(t_nextnb *nnb, char *s)
149 {
150     int i, j, k;
151
152     if (debug)
153     {
154         fprintf(debug, "%s\n", s);
155         fprintf(debug, "nnb->nr: %d\n", nnb->nr);
156         fprintf(debug, "nnb->nrex: %d\n", nnb->nrex);
157         for (i = 0; (i < nnb->nr); i++)
158         {
159             for (j = 0; (j <= nnb->nrex); j++)
160             {
161                 fprintf(debug, "nrexcl[%d][%d]: %d, excl: ", i, j, nnb->nrexcl[i][j]);
162                 for (k = 0; (k < nnb->nrexcl[i][j]); k++)
163                 {
164                     fprintf(debug, "%d, ", nnb->a[i][j][k]);
165                 }
166                 fprintf(debug, "\n");
167             }
168         }
169     }
170 }
171 #endif
172
173 static void nnb2excl(t_nextnb *nnb, t_blocka *excl)
174 {
175     int       i, j, j_index;
176     int       nre, nrx, nrs, nr_of_sortables;
177     sortable *s;
178
179     srenew(excl->index, nnb->nr+1);
180     excl->index[0] = 0;
181     for (i = 0; (i < nnb->nr); i++)
182     {
183         /* calculate the total number of exclusions for atom i */
184         nr_of_sortables = 0;
185         for (nre = 0; (nre <= nnb->nrex); nre++)
186         {
187             nr_of_sortables += nnb->nrexcl[i][nre];
188         }
189
190         if (debug)
191         {
192             fprintf(debug, "nr_of_sortables: %d\n", nr_of_sortables);
193         }
194         /* make space for sortable array */
195         snew(s, nr_of_sortables);
196
197         /* fill the sortable array and sort it */
198         nrs = 0;
199         for (nre = 0; (nre <= nnb->nrex); nre++)
200         {
201             for (nrx = 0; (nrx < nnb->nrexcl[i][nre]); nrx++)
202             {
203                 s[nrs].ai = i;
204                 s[nrs].aj = nnb->a[i][nre][nrx];
205                 nrs++;
206             }
207         }
208         if (nrs != nr_of_sortables)
209         {
210             gmx_incons("Generating exclusions");
211         }
212         prints("nnb2excl before qsort", nr_of_sortables, s);
213         if (nr_of_sortables > 1)
214         {
215             qsort ((void *)s, nr_of_sortables, (size_t)sizeof(s[0]), bond_sort);
216             prints("nnb2excl after qsort", nr_of_sortables, s);
217         }
218
219         /* remove duplicate entries from the list */
220         j_index = 0;
221         if (nr_of_sortables > 0)
222         {
223             for (j = 1; (j < nr_of_sortables); j++)
224             {
225                 if ((s[j].ai != s[j-1].ai) || (s[j].aj != s[j-1].aj))
226                 {
227                     s[j_index++] = s[j-1];
228                 }
229             }
230             s[j_index++] = s[j-1];
231         }
232         nr_of_sortables = j_index;
233         prints("after rm-double", j_index, s);
234
235         /* make space for arrays */
236         srenew(excl->a, excl->nra+nr_of_sortables);
237
238         /* put the sorted exclusions in the target list */
239         for (nrs = 0; (nrs < nr_of_sortables); nrs++)
240         {
241             excl->a[excl->nra+nrs] = s[nrs].aj;
242         }
243         excl->nra       += nr_of_sortables;
244         excl->index[i+1] = excl->nra;
245
246         /* cleanup temporary space */
247         sfree (s);
248     }
249     if (debug)
250     {
251         print_blocka(debug, "Exclusions", "Atom", "Excluded", excl);
252     }
253 }
254
255 static void do_gen(int       nrbonds, /* total number of bonds in s     */
256                    sortable *s,       /* bidirectional list of bonds    */
257                    t_nextnb *nnb)     /* the tmp storage for excl     */
258 /* Assume excl is initalised and s[] contains all bonds bidirectional */
259 {
260     int i, j, k, n, nb;
261
262     /* exclude self */
263     for (i = 0; (i < nnb->nr); i++)
264     {
265         add_nnb(nnb, 0, i, i);
266     }
267     print_nnb(nnb, "After exclude self");
268
269     /* exclude all the bonded atoms */
270     if (nnb->nrex > 0)
271     {
272         for (i = 0; (i < nrbonds); i++)
273         {
274             add_nnb(nnb, 1, s[i].ai, s[i].aj);
275         }
276     }
277     print_nnb(nnb, "After exclude bonds");
278
279     /* for the nr of exclusions per atom */
280     for (n = 1; (n < nnb->nrex); n++)
281     {
282         /* now for all atoms */
283         for (i = 0; (i < nnb->nr); i++)
284         {
285             /* for all directly bonded atoms of atom i */
286             for (j = 0; (j < nnb->nrexcl[i][1]); j++)
287             {
288
289                 /* store the 1st neighbour in nb */
290                 nb = nnb->a[i][1][j];
291
292                 /* store all atoms in nb's n-th list into i's n+1-th list */
293                 for (k = 0; (k < nnb->nrexcl[nb][n]); k++)
294                 {
295                     if (i != nnb->a[nb][n][k])
296                     {
297                         add_nnb(nnb, n+1, i, nnb->a[nb][n][k]);
298                     }
299                 }
300             }
301         }
302     }
303     print_nnb(nnb, "After exclude rest");
304
305 }
306
307 static void add_b(t_params *bonds, int *nrf, sortable *s)
308 {
309     int i;
310     int ai, aj;
311
312     for (i = 0; (i < bonds->nr); i++)
313     {
314         ai = bonds->param[i].AI;
315         aj = bonds->param[i].AJ;
316         if ((ai < 0) || (aj < 0))
317         {
318             gmx_fatal(FARGS, "Impossible atom numbers in bond %d: ai=%d, aj=%d",
319                       i, ai, aj);
320         }
321         /* Add every bond twice */
322         s[(*nrf)].ai   = ai;
323         s[(*nrf)++].aj = aj;
324         s[(*nrf)].aj   = ai;
325         s[(*nrf)++].ai = aj;
326     }
327 }
328
329 void gen_nnb(t_nextnb *nnb, t_params plist[])
330 {
331     sortable *s;
332     int       i, nrbonds, nrf;
333
334     nrbonds = 0;
335     for (i = 0; (i < F_NRE); i++)
336     {
337         if (IS_CHEMBOND(i))
338         {
339             /* we need every bond twice (bidirectional) */
340             nrbonds += 2*plist[i].nr;
341         }
342     }
343
344     snew(s, nrbonds);
345
346     nrf = 0;
347     for (i = 0; (i < F_NRE); i++)
348     {
349         if (IS_CHEMBOND(i))
350         {
351             add_b(&plist[i], &nrf, s);
352         }
353     }
354
355     /* now sort the bonds */
356     prints("gen_excl before qsort", nrbonds, s);
357     if (nrbonds > 1)
358     {
359         qsort((void *) s, nrbonds, (size_t)sizeof(sortable), bond_sort);
360         prints("gen_excl after qsort", nrbonds, s);
361     }
362
363     do_gen(nrbonds, s, nnb);
364     sfree(s);
365 }
366
367 static void
368 sort_and_purge_nnb(t_nextnb *nnb)
369 {
370     int i, j, k, m, n, cnt, found, prev, idx;
371
372     for (i = 0; (i < nnb->nr); i++)
373     {
374         for (n = 0; (n <= nnb->nrex); n++)
375         {
376             /* Sort atoms in this list */
377             qsort(nnb->a[i][n], nnb->nrexcl[i][n], sizeof(int), compare_int);
378
379             cnt  = 0;
380             prev = -1;
381             for (j = 0; j < nnb->nrexcl[i][n]; j++)
382             {
383                 idx = nnb->a[i][n][j];
384
385                 found = 0;
386                 for (m = 0; m < n && !found; m++)
387                 {
388                     for (k = 0; k < nnb->nrexcl[i][m] && !found; k++)
389                     {
390                         found = (idx == nnb->a[i][m][k]);
391                     }
392                 }
393
394                 if (!found && nnb->a[i][n][j] != prev)
395                 {
396                     nnb->a[i][n][cnt] = nnb->a[i][n][j];
397                     prev              = nnb->a[i][n][cnt];
398                     cnt++;
399                 }
400             }
401             nnb->nrexcl[i][n] = cnt;
402         }
403     }
404 }
405
406
407 void generate_excl (int nrexcl, int nratoms, t_params plist[], t_nextnb *nnb, t_blocka *excl)
408 {
409     int i, j, k;
410
411     if (nrexcl < 0)
412     {
413         gmx_fatal(FARGS, "Can't have %d exclusions...", nrexcl);
414     }
415     init_nnb(nnb, nratoms, nrexcl);
416     gen_nnb(nnb, plist);
417     excl->nr = nratoms;
418     sort_and_purge_nnb(nnb);
419     nnb2excl (nnb, excl);
420 }