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