Update copyright statements and change license to LGPL
[alexxy/gromacs.git] / src / kernel / topexcl.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  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 /* This file is completely threadsafe - keep it that way! */
39 #ifdef HAVE_CONFIG_H
40 #include <config.h>
41 #endif
42
43 #include "sysstuff.h"
44 #include "smalloc.h"
45 #include "macros.h"
46 /* #define DEBUG_NNB */
47 #include "gpp_nextnb.h"
48 #include "gmx_fatal.h"
49 #include "toputil.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     return (sa->aj-sb->aj);
65   else
66     return (sa->ai-sb->ai);
67 }
68
69 static int 
70 compare_int (const void * a, const void * b)
71 {
72     return ( *(int*)a - *(int*)b );
73 }
74
75
76 #ifdef DEBUG
77 #define prints(str, n, s) __prints(str, n, s)
78 static void __prints(char *str, int n, sortable *s)
79 {
80   int i;
81
82   if (debug) {
83     fprintf(debug,"%s\n",str);
84     fprintf(debug,"Sortables \n");
85     for (i=0; (i < n); i++)
86       fprintf(debug,"%d\t%d\n",s[i].ai,s[i].aj);
87     
88     fflush(debug);
89   }
90 }
91 #else
92 #define prints(str, n, s)
93 #endif
94
95 void init_nnb(t_nextnb *nnb, int nr, int nrex)
96 {
97   int i;
98
99   /* initiate nnb */
100   nnb->nr   = nr;
101   nnb->nrex = nrex;
102
103   snew(nnb->a,nr);
104   snew(nnb->nrexcl,nr);
105   for (i=0; (i < nr); i++) {
106     snew(nnb->a[i],nrex+1);
107     snew(nnb->nrexcl[i],nrex+1);
108   }
109 }
110
111 static void add_nnb (t_nextnb *nnb, int nre, int i, int j)
112 {
113   srenew(nnb->a[i][nre],nnb->nrexcl[i][nre]+1);
114   nnb->a[i][nre][nnb->nrexcl[i][nre]] = j;
115   nnb->nrexcl[i][nre]++;
116 }
117
118 void done_nnb (t_nextnb *nnb)
119 {
120   int i,nre;
121   
122   for (i=0; (i < nnb->nr); i++) {
123     for (nre=0; (nre <= nnb->nrex); nre++)
124       if (nnb->nrexcl[i][nre] > 0)
125         sfree (nnb->a[i][nre]);
126     sfree (nnb->nrexcl[i]);
127     sfree (nnb->a[i]);
128   }
129   sfree (nnb->a);
130   sfree (nnb->nrexcl);
131   nnb->nr   = 0;
132   nnb->nrex = 0;
133 }
134
135 #ifdef DEBUG_NNB
136 void __print_nnb(t_nextnb *nnb, char *s)
137 {
138   int i,j,k;
139
140   if(debug) {
141     fprintf(debug,"%s\n",s);
142     fprintf(debug,"nnb->nr: %d\n",nnb->nr);
143     fprintf(debug,"nnb->nrex: %d\n",nnb->nrex);
144     for(i=0; (i<nnb->nr); i++)
145       for(j=0; (j<=nnb->nrex); j++) {
146         fprintf(debug,"nrexcl[%d][%d]: %d, excl: ",i,j,nnb->nrexcl[i][j]);
147         for(k=0; (k<nnb->nrexcl[i][j]); k++)
148           fprintf(debug,"%d, ",nnb->a[i][j][k]);
149         fprintf(debug,"\n");
150       }
151   }
152 }
153 #endif
154
155 void nnb2excl (t_nextnb *nnb, t_blocka *excl)
156 {
157   int i,j,j_index;
158   int nre,nrx,nrs,nr_of_sortables;
159   sortable *s;
160
161   srenew(excl->index, nnb->nr+1);
162   excl->index[0]=0;
163   for (i=0; (i < nnb->nr); i++) {
164     /* calculate the total number of exclusions for atom i */
165     nr_of_sortables=0;
166     for (nre=0; (nre<=nnb->nrex); nre++)
167       nr_of_sortables+=nnb->nrexcl[i][nre];
168
169     if (debug)
170       fprintf(debug,"nr_of_sortables: %d\n",nr_of_sortables);
171     /* make space for sortable array */
172     snew(s,nr_of_sortables);
173
174     /* fill the sortable array and sort it */
175     nrs=0;
176     for (nre=0; (nre <= nnb->nrex); nre++)
177       for (nrx=0; (nrx < nnb->nrexcl[i][nre]); nrx++) {
178         s[nrs].ai = i;
179         s[nrs].aj = nnb->a[i][nre][nrx];
180         nrs++;
181       }
182     if (nrs != nr_of_sortables)
183       gmx_incons("Generating exclusions");
184     prints("nnb2excl before qsort",nr_of_sortables,s);
185     if (nr_of_sortables > 1) {
186       qsort ((void *)s,nr_of_sortables,(size_t)sizeof(s[0]),bond_sort);
187       prints("nnb2excl after qsort",nr_of_sortables,s);
188     }
189
190     /* remove duplicate entries from the list */
191     j_index=0;
192     if (nr_of_sortables > 0) {
193       for (j=1; (j < nr_of_sortables); j++)
194         if ((s[j].ai!=s[j-1].ai) || (s[j].aj!=s[j-1].aj))
195           s[j_index++]=s[j-1];
196       s[j_index++]=s[j-1];
197     }
198     nr_of_sortables=j_index;
199     prints("after rm-double",j_index,s);
200
201     /* make space for arrays */
202     srenew(excl->a,excl->nra+nr_of_sortables);
203
204     /* put the sorted exclusions in the target list */
205     for (nrs=0; (nrs<nr_of_sortables); nrs++)
206       excl->a[excl->nra+nrs] = s[nrs].aj;
207     excl->nra+=nr_of_sortables;
208     excl->index[i+1]=excl->nra;
209
210     /* cleanup temporary space */
211     sfree (s);
212   }
213   if (debug)
214     print_blocka(debug,"Exclusions","Atom","Excluded",excl);
215 }
216
217 static void do_gen(int nrbonds,        /* total number of bonds in s    */
218                    sortable *s,        /* bidirectional list of bonds   */
219                    t_nextnb *nnb)      /* the tmp storage for excl     */
220 /* Assume excl is initalised and s[] contains all bonds bidirectional */
221 {
222     int i,j,k,n,nb;
223     
224     /* exclude self */
225     for(i=0; (i<nnb->nr); i++) 
226     {
227         add_nnb(nnb,0,i,i);
228     }
229     print_nnb(nnb,"After exclude self");
230
231     /* exclude all the bonded atoms */
232     if (nnb->nrex > 0)
233     {
234         for (i=0; (i < nrbonds); i++)
235         {
236             add_nnb(nnb,1,s[i].ai,s[i].aj);
237         }
238     }
239     print_nnb(nnb,"After exclude bonds");
240
241   /* for the nr of exclusions per atom */
242     for (n=1; (n < nnb->nrex); n++) 
243     {
244         /* now for all atoms */
245         for (i=0; (i < nnb->nr); i++)
246         {
247             /* for all directly bonded atoms of atom i */
248             for (j=0; (j < nnb->nrexcl[i][1]); j++) 
249             {
250
251                 /* store the 1st neighbour in nb */
252                 nb = nnb->a[i][1][j];
253                 
254                 /* store all atoms in nb's n-th list into i's n+1-th list */
255                 for (k=0; (k < nnb->nrexcl[nb][n]); k++)
256                 {         
257                     if (i != nnb->a[nb][n][k])
258                     {
259                         add_nnb(nnb,n+1,i,nnb->a[nb][n][k]);
260                     }
261                 }
262             }
263         }
264       }
265   print_nnb(nnb,"After exclude rest");
266     
267 }
268
269 static void add_b(t_params *bonds, int *nrf, sortable *s)
270 {
271   int i;
272   int ai,aj;
273   
274   for (i=0; (i < bonds->nr); i++) {
275     ai = bonds->param[i].AI;
276     aj = bonds->param[i].AJ;
277     if ((ai < 0) || (aj < 0)) 
278       gmx_fatal(FARGS,"Impossible atom numbers in bond %d: ai=%d, aj=%d",
279                   i,ai,aj);
280     /* Add every bond twice */
281     s[(*nrf)].ai   = ai;
282     s[(*nrf)++].aj = aj;
283     s[(*nrf)].aj   = ai;
284     s[(*nrf)++].ai = aj;
285   }
286 }
287
288 void gen_nnb(t_nextnb *nnb,t_params plist[])
289 {
290   sortable *s;
291   int      i,nrbonds,nrf;
292
293   nrbonds=0;
294   for(i=0; (i<F_NRE); i++)
295     if (IS_CHEMBOND(i))
296       /* we need every bond twice (bidirectional) */
297       nrbonds += 2*plist[i].nr;
298   
299   snew(s,nrbonds);
300
301   nrf=0;
302   for(i=0; (i<F_NRE); i++)
303     if (IS_CHEMBOND(i))
304       add_b(&plist[i],&nrf,s);
305
306   /* now sort the bonds */
307   prints("gen_excl before qsort",nrbonds,s);
308   if (nrbonds > 1) {
309     qsort((void *) s,nrbonds,(size_t)sizeof(sortable),bond_sort);
310     prints("gen_excl after qsort",nrbonds,s);
311   }
312
313   do_gen(nrbonds,s,nnb);
314   sfree(s);
315 }
316
317 static void
318 sort_and_purge_nnb(t_nextnb *nnb)
319 {
320     int i,j,k,m,n,cnt,found,prev,idx;
321     
322     for (i=0; (i < nnb->nr); i++)
323     {
324         for (n=0; (n <= nnb->nrex); n++) 
325         {
326             /* Sort atoms in this list */
327             qsort(nnb->a[i][n],nnb->nrexcl[i][n],sizeof(int),compare_int);
328             
329             cnt=0;
330             prev=-1;
331             for(j=0;j<nnb->nrexcl[i][n];j++)
332             {
333                 idx = nnb->a[i][n][j];
334                 
335                 found=0;
336                 for(m=0;m<n && !found;m++)
337                 {
338                     for(k=0;k<nnb->nrexcl[i][m] && !found;k++)
339                     {
340                         found = (idx==nnb->a[i][m][k]);
341                     }
342                 }
343                 
344                 if(!found && nnb->a[i][n][j]!=prev)
345                 {
346                     nnb->a[i][n][cnt]=nnb->a[i][n][j];
347                     prev = nnb->a[i][n][cnt];
348                     cnt++;
349                 }
350             }
351             nnb->nrexcl[i][n]=cnt;
352         }
353     }
354 }
355
356
357 void generate_excl (int nrexcl,int nratoms,t_params plist[],t_nextnb *nnb,t_blocka *excl)
358 {
359     int i,j,k;
360     
361     if (nrexcl < 0)
362     {
363         gmx_fatal(FARGS,"Can't have %d exclusions...",nrexcl);
364     }
365     init_nnb(nnb,nratoms,nrexcl);
366     gen_nnb(nnb,plist);
367     excl->nr=nratoms;
368     sort_and_purge_nnb(nnb);
369     nnb2excl (nnb,excl);
370 }
371