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