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