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