41a026d072c18bcb9915c6da9dd162de7cf47132
[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 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <stdlib.h>
43
44 #include "macros.h"
45 /* #define DEBUG_NNB */
46 #include "gpp_nextnb.h"
47 #include "toputil.h"
48
49 #include "gromacs/utility/fatalerror.h"
50 #include "gromacs/utility/smalloc.h"
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         if (debug)
193         {
194             fprintf(debug, "nr_of_sortables: %d\n", nr_of_sortables);
195         }
196         /* make space for sortable array */
197         snew(s, nr_of_sortables);
198
199         /* fill the sortable array and sort it */
200         nrs = 0;
201         for (nre = 0; (nre <= nnb->nrex); nre++)
202         {
203             for (nrx = 0; (nrx < nnb->nrexcl[i][nre]); nrx++)
204             {
205                 s[nrs].ai = i;
206                 s[nrs].aj = nnb->a[i][nre][nrx];
207                 nrs++;
208             }
209         }
210         if (nrs != nr_of_sortables)
211         {
212             gmx_incons("Generating exclusions");
213         }
214         prints("nnb2excl before qsort", nr_of_sortables, s);
215         if (nr_of_sortables > 1)
216         {
217             qsort ((void *)s, nr_of_sortables, (size_t)sizeof(s[0]), bond_sort);
218             prints("nnb2excl after qsort", nr_of_sortables, s);
219         }
220
221         /* remove duplicate entries from the list */
222         j_index = 0;
223         if (nr_of_sortables > 0)
224         {
225             for (j = 1; (j < nr_of_sortables); j++)
226             {
227                 if ((s[j].ai != s[j-1].ai) || (s[j].aj != s[j-1].aj))
228                 {
229                     s[j_index++] = s[j-1];
230                 }
231             }
232             s[j_index++] = s[j-1];
233         }
234         nr_of_sortables = j_index;
235         prints("after rm-double", j_index, s);
236
237         /* make space for arrays */
238         srenew(excl->a, excl->nra+nr_of_sortables);
239
240         /* put the sorted exclusions in the target list */
241         for (nrs = 0; (nrs < nr_of_sortables); nrs++)
242         {
243             excl->a[excl->nra+nrs] = s[nrs].aj;
244         }
245         excl->nra       += nr_of_sortables;
246         excl->index[i+1] = excl->nra;
247
248         /* cleanup temporary space */
249         sfree (s);
250     }
251     if (debug)
252     {
253         print_blocka(debug, "Exclusions", "Atom", "Excluded", excl);
254     }
255 }
256
257 static void do_gen(int       nrbonds, /* total number of bonds in s     */
258                    sortable *s,       /* bidirectional list of bonds    */
259                    t_nextnb *nnb)     /* the tmp storage for excl     */
260 /* Assume excl is initalised and s[] contains all bonds bidirectional */
261 {
262     int i, j, k, n, nb;
263
264     /* exclude self */
265     for (i = 0; (i < nnb->nr); i++)
266     {
267         add_nnb(nnb, 0, i, i);
268     }
269     print_nnb(nnb, "After exclude self");
270
271     /* exclude all the bonded atoms */
272     if (nnb->nrex > 0)
273     {
274         for (i = 0; (i < nrbonds); i++)
275         {
276             add_nnb(nnb, 1, s[i].ai, s[i].aj);
277         }
278     }
279     print_nnb(nnb, "After exclude bonds");
280
281     /* for the nr of exclusions per atom */
282     for (n = 1; (n < nnb->nrex); n++)
283     {
284         /* now for all atoms */
285         for (i = 0; (i < nnb->nr); i++)
286         {
287             /* for all directly bonded atoms of atom i */
288             for (j = 0; (j < nnb->nrexcl[i][1]); j++)
289             {
290
291                 /* store the 1st neighbour in nb */
292                 nb = nnb->a[i][1][j];
293
294                 /* store all atoms in nb's n-th list into i's n+1-th list */
295                 for (k = 0; (k < nnb->nrexcl[nb][n]); k++)
296                 {
297                     if (i != nnb->a[nb][n][k])
298                     {
299                         add_nnb(nnb, n+1, i, nnb->a[nb][n][k]);
300                     }
301                 }
302             }
303         }
304     }
305     print_nnb(nnb, "After exclude rest");
306
307 }
308
309 static void add_b(t_params *bonds, int *nrf, sortable *s)
310 {
311     int i;
312     int ai, aj;
313
314     for (i = 0; (i < bonds->nr); i++)
315     {
316         ai = bonds->param[i].AI;
317         aj = bonds->param[i].AJ;
318         if ((ai < 0) || (aj < 0))
319         {
320             gmx_fatal(FARGS, "Impossible atom numbers in bond %d: ai=%d, aj=%d",
321                       i, ai, aj);
322         }
323         /* Add every bond twice */
324         s[(*nrf)].ai   = ai;
325         s[(*nrf)++].aj = aj;
326         s[(*nrf)].aj   = ai;
327         s[(*nrf)++].ai = aj;
328     }
329 }
330
331 void gen_nnb(t_nextnb *nnb, t_params plist[])
332 {
333     sortable *s;
334     int       i, nrbonds, nrf;
335
336     nrbonds = 0;
337     for (i = 0; (i < F_NRE); i++)
338     {
339         if (IS_CHEMBOND(i))
340         {
341             /* we need every bond twice (bidirectional) */
342             nrbonds += 2*plist[i].nr;
343         }
344     }
345
346     snew(s, nrbonds);
347
348     nrf = 0;
349     for (i = 0; (i < F_NRE); i++)
350     {
351         if (IS_CHEMBOND(i))
352         {
353             add_b(&plist[i], &nrf, s);
354         }
355     }
356
357     /* now sort the bonds */
358     prints("gen_excl before qsort", nrbonds, s);
359     if (nrbonds > 1)
360     {
361         qsort((void *) s, nrbonds, (size_t)sizeof(sortable), bond_sort);
362         prints("gen_excl after qsort", nrbonds, s);
363     }
364
365     do_gen(nrbonds, s, nnb);
366     sfree(s);
367 }
368
369 static void
370 sort_and_purge_nnb(t_nextnb *nnb)
371 {
372     int i, j, k, m, n, cnt, found, prev, idx;
373
374     for (i = 0; (i < nnb->nr); i++)
375     {
376         for (n = 0; (n <= nnb->nrex); n++)
377         {
378             /* Sort atoms in this list */
379             qsort(nnb->a[i][n], nnb->nrexcl[i][n], sizeof(int), compare_int);
380
381             cnt  = 0;
382             prev = -1;
383             for (j = 0; j < nnb->nrexcl[i][n]; j++)
384             {
385                 idx = nnb->a[i][n][j];
386
387                 found = 0;
388                 for (m = 0; m < n && !found; m++)
389                 {
390                     for (k = 0; k < nnb->nrexcl[i][m] && !found; k++)
391                     {
392                         found = (idx == nnb->a[i][m][k]);
393                     }
394                 }
395
396                 if (!found && nnb->a[i][n][j] != prev)
397                 {
398                     nnb->a[i][n][cnt] = nnb->a[i][n][j];
399                     prev              = nnb->a[i][n][cnt];
400                     cnt++;
401                 }
402             }
403             nnb->nrexcl[i][n] = cnt;
404         }
405     }
406 }
407
408
409 void generate_excl (int nrexcl, int nratoms, t_params plist[], t_nextnb *nnb, t_blocka *excl)
410 {
411     int i, j, k;
412
413     if (nrexcl < 0)
414     {
415         gmx_fatal(FARGS, "Can't have %d exclusions...", nrexcl);
416     }
417     init_nnb(nnb, nratoms, nrexcl);
418     gen_nnb(nnb, plist);
419     excl->nr = nratoms;
420     sort_and_purge_nnb(nnb);
421     nnb2excl (nnb, excl);
422 }