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