Sort all includes in src/gromacs
[alexxy/gromacs.git] / src / gromacs / gmxlib / splitter.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) 2013,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 #include "gmxpre.h"
38
39 #include "gromacs/legacyheaders/splitter.h"
40
41 #include <assert.h>
42 #include <stdlib.h>
43 #include <string.h>
44
45 #include "gromacs/legacyheaders/macros.h"
46 #include "gromacs/pbcutil/mshift.h"
47 #include "gromacs/topology/block.h"
48 #include "gromacs/topology/idef.h"
49 #include "gromacs/utility/fatalerror.h"
50 #include "gromacs/utility/smalloc.h"
51
52 typedef struct {
53     int atom, sid;
54 } t_sid;
55
56 static int sid_comp(const void *a, const void *b)
57 {
58     t_sid *sa, *sb;
59     int    dd;
60
61     sa = (t_sid *)a;
62     sb = (t_sid *)b;
63
64     dd = sa->sid-sb->sid;
65     if (dd == 0)
66     {
67         return (sa->atom-sb->atom);
68     }
69     else
70     {
71         return dd;
72     }
73 }
74
75 static int mk_grey(egCol egc[], t_graph *g, int *AtomI,
76                    int maxsid, t_sid sid[])
77 {
78     int  j, ng, ai, aj, g0;
79
80     ng = 0;
81     ai = *AtomI;
82
83     g0 = g->at_start;
84     /* Loop over all the bonds */
85     for (j = 0; (j < g->nedge[ai]); j++)
86     {
87         aj = g->edge[ai][j]-g0;
88         /* If there is a white one, make it gray and set pbc */
89         if (egc[aj] == egcolWhite)
90         {
91             if (aj < *AtomI)
92             {
93                 *AtomI = aj;
94             }
95             egc[aj] = egcolGrey;
96
97             /* Check whether this one has been set before... */
98             range_check(aj+g0, 0, maxsid);
99             range_check(ai+g0, 0, maxsid);
100             if (sid[aj+g0].sid != -1)
101             {
102                 gmx_fatal(FARGS, "sid[%d]=%d, sid[%d]=%d, file %s, line %d",
103                           ai, sid[ai+g0].sid, aj, sid[aj+g0].sid, __FILE__, __LINE__);
104             }
105             else
106             {
107                 sid[aj+g0].sid  = sid[ai+g0].sid;
108                 sid[aj+g0].atom = aj+g0;
109             }
110             ng++;
111         }
112     }
113     return ng;
114 }
115
116 static int first_colour(int fC, egCol Col, t_graph *g, egCol egc[])
117 /* Return the first node with colour Col starting at fC.
118  * return -1 if none found.
119  */
120 {
121     int i;
122
123     for (i = fC; (i < g->nnodes); i++)
124     {
125         if ((g->nedge[i] > 0) && (egc[i] == Col))
126         {
127             return i;
128         }
129     }
130
131     return -1;
132 }
133
134 static int mk_sblocks(FILE *fp, t_graph *g, int maxsid, t_sid sid[])
135 {
136     int     ng, nnodes;
137     int     nW, nG, nB; /* Number of Grey, Black, White */
138     int     fW, fG;     /* First of each category       */
139     egCol  *egc = NULL; /* The colour of each node      */
140     int     g0, nblock;
141
142     if (!g->nbound)
143     {
144         return 0;
145     }
146
147     nblock = 0;
148
149     nnodes = g->nnodes;
150     snew(egc, nnodes);
151
152     g0 = g->at_start;
153     nW = g->nbound;
154     nG = 0;
155     nB = 0;
156
157     fW = 0;
158
159     /* We even have a loop invariant:
160      * nW+nG+nB == g->nbound
161      */
162
163     if (fp)
164     {
165         fprintf(fp, "Walking down the molecule graph to make constraint-blocks\n");
166     }
167
168     while (nW > 0)
169     {
170         /* Find the first white, this will allways be a larger
171          * number than before, because no nodes are made white
172          * in the loop
173          */
174         if ((fW = first_colour(fW, egcolWhite, g, egc)) == -1)
175         {
176             gmx_fatal(FARGS, "No WHITE nodes found while nW=%d\n", nW);
177         }
178
179         /* Make the first white node grey, and set the block number */
180         egc[fW]        = egcolGrey;
181         range_check(fW+g0, 0, maxsid);
182         sid[fW+g0].sid = nblock++;
183         nG++;
184         nW--;
185
186         /* Initial value for the first grey */
187         fG = fW;
188
189         if (debug)
190         {
191             fprintf(debug, "Starting G loop (nW=%d, nG=%d, nB=%d, total %d)\n",
192                     nW, nG, nB, nW+nG+nB);
193         }
194
195         while (nG > 0)
196         {
197             if ((fG = first_colour(fG, egcolGrey, g, egc)) == -1)
198             {
199                 gmx_fatal(FARGS, "No GREY nodes found while nG=%d\n", nG);
200             }
201
202             /* Make the first grey node black */
203             egc[fG] = egcolBlack;
204             nB++;
205             nG--;
206
207             /* Make all the neighbours of this black node grey
208              * and set their block number
209              */
210             ng = mk_grey(egc, g, &fG, maxsid, sid);
211             /* ng is the number of white nodes made grey */
212             nG += ng;
213             nW -= ng;
214         }
215     }
216     sfree(egc);
217
218     if (debug)
219     {
220         fprintf(debug, "Found %d shake blocks\n", nblock);
221     }
222
223     return nblock;
224 }
225
226
227 typedef struct {
228     int first, last, sid;
229 } t_merge_sid;
230
231 static int ms_comp(const void *a, const void *b)
232 {
233     t_merge_sid *ma = (t_merge_sid *)a;
234     t_merge_sid *mb = (t_merge_sid *)b;
235     int          d;
236
237     d = ma->first-mb->first;
238     if (d == 0)
239     {
240         return ma->last-mb->last;
241     }
242     else
243     {
244         return d;
245     }
246 }
247
248 static int merge_sid(int at_start, int at_end, int nsid, t_sid sid[],
249                      t_blocka *sblock)
250 {
251     int          i, j, k, n, isid, ndel;
252     t_merge_sid *ms;
253     int          nChanged;
254
255     /* We try to remdy the following problem:
256      * Atom: 1  2  3  4  5 6 7 8 9 10
257      * Sid:  0 -1  1  0 -1 1 1 1 1 1
258      */
259
260     /* Determine first and last atom in each shake ID */
261     snew(ms, nsid);
262
263     for (k = 0; (k < nsid); k++)
264     {
265         ms[k].first = at_end+1;
266         ms[k].last  = -1;
267         ms[k].sid   = k;
268     }
269     for (i = at_start; (i < at_end); i++)
270     {
271         isid = sid[i].sid;
272         range_check(isid, -1, nsid);
273         if (isid >= 0)
274         {
275             ms[isid].first = min(ms[isid].first, sid[i].atom);
276             ms[isid].last  = max(ms[isid].last, sid[i].atom);
277         }
278     }
279     qsort(ms, nsid, sizeof(ms[0]), ms_comp);
280
281     /* Now merge the overlapping ones */
282     ndel = 0;
283     for (k = 0; (k < nsid); )
284     {
285         for (j = k+1; (j < nsid); )
286         {
287             if (ms[j].first <= ms[k].last)
288             {
289                 ms[k].last  = max(ms[k].last, ms[j].last);
290                 ms[k].first = min(ms[k].first, ms[j].first);
291                 ms[j].sid   = -1;
292                 ndel++;
293                 j++;
294             }
295             else
296             {
297                 k = j;
298                 j = k+1;
299             }
300         }
301         if (j == nsid)
302         {
303             k++;
304         }
305     }
306     for (k = 0; (k < nsid); k++)
307     {
308         while ((k < nsid-1) && (ms[k].sid == -1))
309         {
310             for (j = k+1; (j < nsid); j++)
311             {
312                 memcpy(&(ms[j-1]), &(ms[j]), sizeof(ms[0]));
313             }
314             nsid--;
315         }
316     }
317
318     for (k = at_start; (k < at_end); k++)
319     {
320         sid[k].atom = k;
321         sid[k].sid  = -1;
322     }
323     sblock->nr           = nsid;
324     sblock->nalloc_index = sblock->nr+1;
325     snew(sblock->index, sblock->nalloc_index);
326     sblock->nra      = at_end - at_start;
327     sblock->nalloc_a = sblock->nra;
328     snew(sblock->a, sblock->nalloc_a);
329     sblock->index[0] = 0;
330     for (k = n = 0; (k < nsid); k++)
331     {
332         sblock->index[k+1] = sblock->index[k] + ms[k].last - ms[k].first+1;
333         for (j = ms[k].first; (j <= ms[k].last); j++)
334         {
335             range_check(n, 0, sblock->nra);
336             sblock->a[n++] = j;
337             range_check(j, 0, at_end);
338             if (sid[j].sid == -1)
339             {
340                 sid[j].sid = k;
341             }
342             else
343             {
344                 fprintf(stderr, "Double sids (%d, %d) for atom %d\n", sid[j].sid, k, j);
345             }
346         }
347     }
348     assert(k == nsid);
349     /* Removed 2007-09-04
350        sblock->index[k+1] = natoms;
351        for(k=0; (k<natoms); k++)
352        if (sid[k].sid == -1)
353         sblock->a[n++] = k;
354        assert(n == natoms);
355      */
356     sblock->nra = n;
357     assert(sblock->index[k] == sblock->nra);
358     sfree(ms);
359
360     return nsid;
361 }
362
363 void gen_sblocks(FILE *fp, int at_start, int at_end,
364                  t_idef *idef, t_blocka *sblock,
365                  gmx_bool bSettle)
366 {
367     t_graph *g;
368     int      i, i0, j, k, istart, n;
369     t_sid   *sid;
370     int      isid, nsid;
371
372     g = mk_graph(NULL, idef, at_start, at_end, TRUE, bSettle);
373     if (debug)
374     {
375         p_graph(debug, "Graaf Dracula", g);
376     }
377     snew(sid, at_end);
378     for (i = at_start; (i < at_end); i++)
379     {
380         sid[i].atom =  i;
381         sid[i].sid  = -1;
382     }
383     nsid = mk_sblocks(fp, g, at_end, sid);
384
385     if (!nsid)
386     {
387         return;
388     }
389
390     /* Now sort the shake blocks... */
391     qsort(sid+at_start, at_end-at_start, (size_t)sizeof(sid[0]), sid_comp);
392
393     if (debug)
394     {
395         fprintf(debug, "Sorted shake block\n");
396         for (i = at_start; (i < at_end); i++)
397         {
398             fprintf(debug, "sid[%5d] = atom:%5d sid:%5d\n", i, sid[i].atom, sid[i].sid);
399         }
400     }
401     /* Now check how many are NOT -1, i.e. how many have to be shaken */
402     for (i0 = at_start; (i0 < at_end); i0++)
403     {
404         if (sid[i0].sid > -1)
405         {
406             break;
407         }
408     }
409
410     /* Now we have the sids that have to be shaken. We'll check the min and
411      * max atom numbers and this determines the shake block. DvdS 2007-07-19.
412      * For the purpose of making boundaries all atoms in between need to be
413      * part of the shake block too. There may be cases where blocks overlap
414      * and they will have to be merged.
415      */
416     nsid = merge_sid(at_start, at_end, nsid, sid, sblock);
417     /* Now sort the shake blocks again... */
418     /*qsort(sid,natoms,(size_t)sizeof(sid[0]),sid_comp);*/
419
420     /* Fill the sblock struct */
421     /*  sblock->nr  = nsid;
422        sblock->nra = natoms;
423        srenew(sblock->a,sblock->nra);
424        srenew(sblock->index,sblock->nr+1);
425
426        i    = i0;
427        isid = sid[i].sid;
428        n    = k = 0;
429        sblock->index[n++]=k;
430        while (i < natoms) {
431        istart = sid[i].atom;
432        while ((i<natoms-1) && (sid[i+1].sid == isid))
433        i++;*/
434     /* After while: we found a new block, or are thru with the atoms */
435     /*    for(j=istart; (j<=sid[i].atom); j++,k++)
436         sblock->a[k]=j;
437        sblock->index[n] = k;
438        if (i < natoms-1)
439         n++;
440        if (n > nsid)
441         gmx_fatal(FARGS,"Death Horror: nsid = %d, n= %d",nsid,n);
442        i++;
443        isid = sid[i].sid;
444        }
445      */
446     sfree(sid);
447     /* Due to unknown reason this free generates a problem sometimes */
448     done_graph(g);
449     sfree(g);
450     if (debug)
451     {
452         fprintf(debug, "Done gen_sblocks\n");
453     }
454 }