1bbbd9ed097093364b27d4a5f93f03dad2dea898
[alexxy/gromacs.git] / src / gromacs / mdlib / splitter.cpp
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,2015, 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 "splitter.h"
40
41 #include <cstdlib>
42 #include <cstring>
43
44 #include <algorithm>
45
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/gmxassert.h"
51 #include "gromacs/utility/smalloc.h"
52
53 typedef struct {
54     int atom, sid;
55 } t_sid;
56
57 static int sid_comp(const void *a, const void *b)
58 {
59     t_sid *sa, *sb;
60     int    dd;
61
62     sa = (t_sid *)a;
63     sb = (t_sid *)b;
64
65     dd = sa->sid-sb->sid;
66     if (dd == 0)
67     {
68         return (sa->atom-sb->atom);
69     }
70     else
71     {
72         return dd;
73     }
74 }
75
76 static int mk_grey(egCol egc[], t_graph *g, int *AtomI,
77                    int maxsid, t_sid sid[])
78 {
79     int  j, ng, ai, aj, g0;
80
81     ng = 0;
82     ai = *AtomI;
83
84     g0 = g->at_start;
85     /* Loop over all the bonds */
86     for (j = 0; (j < g->nedge[ai]); j++)
87     {
88         aj = g->edge[ai][j]-g0;
89         /* If there is a white one, make it gray and set pbc */
90         if (egc[aj] == egcolWhite)
91         {
92             if (aj < *AtomI)
93             {
94                 *AtomI = aj;
95             }
96             egc[aj] = egcolGrey;
97
98             /* Check whether this one has been set before... */
99             range_check(aj+g0, 0, maxsid);
100             range_check(ai+g0, 0, maxsid);
101             if (sid[aj+g0].sid != -1)
102             {
103                 gmx_fatal(FARGS, "sid[%d]=%d, sid[%d]=%d, file %s, line %d",
104                           ai, sid[ai+g0].sid, aj, sid[aj+g0].sid, __FILE__, __LINE__);
105             }
106             else
107             {
108                 sid[aj+g0].sid  = sid[ai+g0].sid;
109                 sid[aj+g0].atom = aj+g0;
110             }
111             ng++;
112         }
113     }
114     return ng;
115 }
116
117 static int first_colour(int fC, egCol Col, t_graph *g, egCol egc[])
118 /* Return the first node with colour Col starting at fC.
119  * return -1 if none found.
120  */
121 {
122     int i;
123
124     for (i = fC; (i < g->nnodes); i++)
125     {
126         if ((g->nedge[i] > 0) && (egc[i] == Col))
127         {
128             return i;
129         }
130     }
131
132     return -1;
133 }
134
135 static int mk_sblocks(FILE *fp, t_graph *g, int maxsid, t_sid sid[])
136 {
137     int     ng, nnodes;
138     int     nW, nG, nB; /* Number of Grey, Black, White */
139     int     fW, fG;     /* First of each category       */
140     egCol  *egc = NULL; /* The colour of each node      */
141     int     g0, nblock;
142
143     if (!g->nbound)
144     {
145         return 0;
146     }
147
148     nblock = 0;
149
150     nnodes = g->nnodes;
151     snew(egc, nnodes);
152
153     g0 = g->at_start;
154     nW = g->nbound;
155     nG = 0;
156     nB = 0;
157
158     fW = 0;
159
160     /* We even have a loop invariant:
161      * nW+nG+nB == g->nbound
162      */
163
164     if (fp)
165     {
166         fprintf(fp, "Walking down the molecule graph to make constraint-blocks\n");
167     }
168
169     while (nW > 0)
170     {
171         /* Find the first white, this will allways be a larger
172          * number than before, because no nodes are made white
173          * in the loop
174          */
175         if ((fW = first_colour(fW, egcolWhite, g, egc)) == -1)
176         {
177             gmx_fatal(FARGS, "No WHITE nodes found while nW=%d\n", nW);
178         }
179
180         /* Make the first white node grey, and set the block number */
181         egc[fW]        = egcolGrey;
182         range_check(fW+g0, 0, maxsid);
183         sid[fW+g0].sid = nblock++;
184         nG++;
185         nW--;
186
187         /* Initial value for the first grey */
188         fG = fW;
189
190         if (debug)
191         {
192             fprintf(debug, "Starting G loop (nW=%d, nG=%d, nB=%d, total %d)\n",
193                     nW, nG, nB, nW+nG+nB);
194         }
195
196         while (nG > 0)
197         {
198             if ((fG = first_colour(fG, egcolGrey, g, egc)) == -1)
199             {
200                 gmx_fatal(FARGS, "No GREY nodes found while nG=%d\n", nG);
201             }
202
203             /* Make the first grey node black */
204             egc[fG] = egcolBlack;
205             nB++;
206             nG--;
207
208             /* Make all the neighbours of this black node grey
209              * and set their block number
210              */
211             ng = mk_grey(egc, g, &fG, maxsid, sid);
212             /* ng is the number of white nodes made grey */
213             nG += ng;
214             nW -= ng;
215         }
216     }
217     sfree(egc);
218
219     if (debug)
220     {
221         fprintf(debug, "Found %d shake blocks\n", nblock);
222     }
223
224     return nblock;
225 }
226
227
228 typedef struct {
229     int first, last, sid;
230 } t_merge_sid;
231
232 static int ms_comp(const void *a, const void *b)
233 {
234     const t_merge_sid *ma = reinterpret_cast<const t_merge_sid *>(a);
235     const t_merge_sid *mb = reinterpret_cast<const t_merge_sid *>(b);
236     int                d;
237
238     d = ma->first-mb->first;
239     if (d == 0)
240     {
241         return ma->last-mb->last;
242     }
243     else
244     {
245         return d;
246     }
247 }
248
249 static int merge_sid(int at_start, int at_end, int nsid, t_sid sid[],
250                      t_blocka *sblock)
251 {
252     int          i, j, k, n, isid, ndel;
253     t_merge_sid *ms;
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 = std::min(ms[isid].first, sid[i].atom);
276             ms[isid].last  = std::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  = std::max(ms[k].last, ms[j].last);
290                 ms[k].first = std::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                 std::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
349     sblock->nra = n;
350     GMX_RELEASE_ASSERT(sblock->index[k] == sblock->nra, "Internal inconsistency; sid merge failed");
351
352     sfree(ms);
353
354     return nsid;
355 }
356
357 void gen_sblocks(FILE *fp, int at_start, int at_end,
358                  t_idef *idef, t_blocka *sblock,
359                  gmx_bool bSettle)
360 {
361     t_graph *g;
362     int      i, i0;
363     t_sid   *sid;
364     int      nsid;
365
366     g = mk_graph(NULL, idef, at_start, at_end, TRUE, bSettle);
367     if (debug)
368     {
369         p_graph(debug, "Graaf Dracula", g);
370     }
371     snew(sid, at_end);
372     for (i = at_start; (i < at_end); i++)
373     {
374         sid[i].atom =  i;
375         sid[i].sid  = -1;
376     }
377     nsid = mk_sblocks(fp, g, at_end, sid);
378
379     if (!nsid)
380     {
381         return;
382     }
383
384     /* Now sort the shake blocks... */
385     qsort(sid+at_start, at_end-at_start, static_cast<size_t>(sizeof(sid[0])), sid_comp);
386
387     if (debug)
388     {
389         fprintf(debug, "Sorted shake block\n");
390         for (i = at_start; (i < at_end); i++)
391         {
392             fprintf(debug, "sid[%5d] = atom:%5d sid:%5d\n", i, sid[i].atom, sid[i].sid);
393         }
394     }
395     /* Now check how many are NOT -1, i.e. how many have to be shaken */
396     for (i0 = at_start; (i0 < at_end); i0++)
397     {
398         if (sid[i0].sid > -1)
399         {
400             break;
401         }
402     }
403
404     /* Now we have the sids that have to be shaken. We'll check the min and
405      * max atom numbers and this determines the shake block. DvdS 2007-07-19.
406      * For the purpose of making boundaries all atoms in between need to be
407      * part of the shake block too. There may be cases where blocks overlap
408      * and they will have to be merged.
409      */
410     merge_sid(at_start, at_end, nsid, sid, sblock);
411     sfree(sid);
412     /* Due to unknown reason this free generates a problem sometimes */
413     done_graph(g);
414     sfree(g);
415     if (debug)
416     {
417         fprintf(debug, "Done gen_sblocks\n");
418     }
419 }