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