a970c3fc0a71a7ba54ef9d67deb8d00f391f9cb2
[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 by the GROMACS development team.
7  * Copyright (c) 2019,2020, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #include "gmxpre.h"
39
40 #include "splitter.h"
41
42 #include <cstdlib>
43 #include <cstring>
44
45 #include <algorithm>
46
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 {
56     int atom, sid;
57 } t_sid;
58
59 static bool sid_comp(const t_sid& sa, const t_sid& sb)
60 {
61     if (sa.sid == sb.sid)
62     {
63         return sa.atom < sb.atom;
64     }
65     else
66     {
67         return sa.sid < sb.sid;
68     }
69 }
70
71 static int mk_grey(gmx::ArrayRef<egCol> edgeColor, const t_graph* g, int* AtomI, int maxsid, t_sid sid[])
72 {
73     int ng, ai, g0;
74
75     ng = 0;
76     ai = *AtomI;
77
78     g0 = g->edgeAtomBegin;
79     /* Loop over all the bonds */
80     for (int aj : g->edges[ai])
81     {
82         aj -= g0;
83         /* If there is a white one, make it gray and set pbc */
84         if (edgeColor[aj] == egcolWhite)
85         {
86             if (aj < *AtomI)
87             {
88                 *AtomI = aj;
89             }
90             edgeColor[aj] = egcolGrey;
91
92             /* Check whether this one has been set before... */
93             range_check(aj + g0, 0, maxsid);
94             range_check(ai + g0, 0, maxsid);
95             if (sid[aj + g0].sid != -1)
96             {
97                 gmx_fatal(FARGS, "sid[%d]=%d, sid[%d]=%d, file %s, line %d", ai, sid[ai + g0].sid,
98                           aj, sid[aj + g0].sid, __FILE__, __LINE__);
99             }
100             else
101             {
102                 sid[aj + g0].sid  = sid[ai + g0].sid;
103                 sid[aj + g0].atom = aj + g0;
104             }
105             ng++;
106         }
107     }
108     return ng;
109 }
110
111 static int first_colour(const int fC, const egCol Col, const t_graph* g, gmx::ArrayRef<const egCol> edgeColor)
112 /* Return the first node with colour Col starting at fC.
113  * return -1 if none found.
114  */
115 {
116     int i;
117
118     for (i = fC; i < int(g->edges.size()); i++)
119     {
120         if (!g->edges[i].empty() && edgeColor[i] == Col)
121         {
122             return i;
123         }
124     }
125
126     return -1;
127 }
128
129 static int mk_sblocks(FILE* fp, t_graph* g, int maxsid, t_sid sid[])
130 {
131     int ng;
132     int nW, nG, nB; /* Number of Grey, Black, White     */
133     int fW, fG;     /* First of each category   */
134     int g0, nblock;
135
136     if (!g->numConnectedAtoms)
137     {
138         return 0;
139     }
140
141     nblock = 0;
142
143     std::vector<egCol> edgeColor(g->edges.size(), egcolWhite);
144
145     g0 = g->edgeAtomBegin;
146     nW = g->numConnectedAtoms;
147     nG = 0;
148     nB = 0;
149
150     fW = 0;
151
152     /* We even have a loop invariant:
153      * nW+nG+nB == g->nbound
154      */
155
156     if (fp)
157     {
158         fprintf(fp, "Walking down the molecule graph to make constraint-blocks\n");
159     }
160
161     while (nW > 0)
162     {
163         /* Find the first white, this will allways be a larger
164          * number than before, because no nodes are made white
165          * in the loop
166          */
167         if ((fW = first_colour(fW, egcolWhite, g, edgeColor)) == -1)
168         {
169             gmx_fatal(FARGS, "No WHITE nodes found while nW=%d\n", nW);
170         }
171
172         /* Make the first white node grey, and set the block number */
173         edgeColor[fW] = egcolGrey;
174         range_check(fW + g0, 0, maxsid);
175         sid[fW + g0].sid = nblock++;
176         nG++;
177         nW--;
178
179         /* Initial value for the first grey */
180         fG = fW;
181
182         if (debug)
183         {
184             fprintf(debug, "Starting G loop (nW=%d, nG=%d, nB=%d, total %d)\n", nW, nG, nB, nW + nG + nB);
185         }
186
187         while (nG > 0)
188         {
189             if ((fG = first_colour(fG, egcolGrey, g, edgeColor)) == -1)
190             {
191                 gmx_fatal(FARGS, "No GREY nodes found while nG=%d\n", nG);
192             }
193
194             /* Make the first grey node black */
195             edgeColor[fG] = egcolBlack;
196             nB++;
197             nG--;
198
199             /* Make all the neighbours of this black node grey
200              * and set their block number
201              */
202             ng = mk_grey(edgeColor, g, &fG, maxsid, sid);
203             /* ng is the number of white nodes made grey */
204             nG += ng;
205             nW -= ng;
206         }
207     }
208
209     if (debug)
210     {
211         fprintf(debug, "Found %d shake blocks\n", nblock);
212     }
213
214     return nblock;
215 }
216
217
218 typedef struct
219 {
220     int first, last, sid;
221 } t_merge_sid;
222
223 static int ms_comp(const void* a, const void* b)
224 {
225     const t_merge_sid* ma = reinterpret_cast<const t_merge_sid*>(a);
226     const t_merge_sid* mb = reinterpret_cast<const t_merge_sid*>(b);
227     int                d;
228
229     d = ma->first - mb->first;
230     if (d == 0)
231     {
232         return ma->last - mb->last;
233     }
234     else
235     {
236         return d;
237     }
238 }
239
240 static int merge_sid(int at_start, int at_end, int nsid, t_sid sid[], t_blocka* sblock)
241 {
242     int          i, j, k, n, isid, ndel;
243     t_merge_sid* ms;
244
245     /* We try to remdy the following problem:
246      * Atom: 1  2  3  4  5 6 7 8 9 10
247      * Sid:  0 -1  1  0 -1 1 1 1 1 1
248      */
249
250     /* Determine first and last atom in each shake ID */
251     snew(ms, nsid);
252
253     for (k = 0; (k < nsid); k++)
254     {
255         ms[k].first = at_end + 1;
256         ms[k].last  = -1;
257         ms[k].sid   = k;
258     }
259     for (i = at_start; (i < at_end); i++)
260     {
261         isid = sid[i].sid;
262         range_check(isid, -1, nsid);
263         if (isid >= 0)
264         {
265             ms[isid].first = std::min(ms[isid].first, sid[i].atom);
266             ms[isid].last  = std::max(ms[isid].last, sid[i].atom);
267         }
268     }
269     qsort(ms, nsid, sizeof(ms[0]), ms_comp);
270
271     /* Now merge the overlapping ones */
272     ndel = 0;
273     for (k = 0; (k < nsid);)
274     {
275         for (j = k + 1; (j < nsid);)
276         {
277             if (ms[j].first <= ms[k].last)
278             {
279                 ms[k].last  = std::max(ms[k].last, ms[j].last);
280                 ms[k].first = std::min(ms[k].first, ms[j].first);
281                 ms[j].sid   = -1;
282                 ndel++;
283                 j++;
284             }
285             else
286             {
287                 k = j;
288                 j = k + 1;
289             }
290         }
291         if (j == nsid)
292         {
293             k++;
294         }
295     }
296     for (k = 0; (k < nsid); k++)
297     {
298         while ((k < nsid - 1) && (ms[k].sid == -1))
299         {
300             for (j = k + 1; (j < nsid); j++)
301             {
302                 std::memcpy(&(ms[j - 1]), &(ms[j]), sizeof(ms[0]));
303             }
304             nsid--;
305         }
306     }
307
308     for (k = at_start; (k < at_end); k++)
309     {
310         sid[k].atom = k;
311         sid[k].sid  = -1;
312     }
313     sblock->nr           = nsid;
314     sblock->nalloc_index = sblock->nr + 1;
315     snew(sblock->index, sblock->nalloc_index);
316     sblock->nra      = at_end - at_start;
317     sblock->nalloc_a = sblock->nra;
318     snew(sblock->a, sblock->nalloc_a);
319     sblock->index[0] = 0;
320     for (k = n = 0; (k < nsid); k++)
321     {
322         sblock->index[k + 1] = sblock->index[k] + ms[k].last - ms[k].first + 1;
323         for (j = ms[k].first; (j <= ms[k].last); j++)
324         {
325             range_check(n, 0, sblock->nra);
326             sblock->a[n++] = j;
327             range_check(j, 0, at_end);
328             if (sid[j].sid == -1)
329             {
330                 sid[j].sid = k;
331             }
332             else
333             {
334                 fprintf(stderr, "Double sids (%d, %d) for atom %d\n", sid[j].sid, k, j);
335             }
336         }
337     }
338
339     sblock->nra = n;
340     GMX_RELEASE_ASSERT(sblock->index[k] == sblock->nra, "Internal inconsistency; sid merge failed");
341
342     sfree(ms);
343
344     return nsid;
345 }
346
347 void gen_sblocks(FILE* fp, int at_end, const InteractionDefinitions& idef, t_blocka* sblock, gmx_bool bSettle)
348 {
349     t_graph* g;
350     int      i, i0;
351     t_sid*   sid;
352     int      nsid;
353
354     g = mk_graph(nullptr, idef, at_end, TRUE, bSettle);
355     if (debug)
356     {
357         p_graph(debug, "Graaf Dracula", g);
358     }
359     snew(sid, at_end);
360     for (i = 0; i < at_end; i++)
361     {
362         sid[i].atom = i;
363         sid[i].sid  = -1;
364     }
365     nsid = mk_sblocks(fp, g, at_end, sid);
366
367     if (!nsid)
368     {
369         return;
370     }
371
372     /* Now sort the shake blocks... */
373     std::sort(sid, sid + at_end, sid_comp);
374
375     if (debug)
376     {
377         fprintf(debug, "Sorted shake block\n");
378         for (i = 0; i < at_end; i++)
379         {
380             fprintf(debug, "sid[%5d] = atom:%5d sid:%5d\n", i, sid[i].atom, sid[i].sid);
381         }
382     }
383     /* Now check how many are NOT -1, i.e. how many have to be shaken */
384     for (i0 = 0; i0 < at_end; i0++)
385     {
386         if (sid[i0].sid > -1)
387         {
388             break;
389         }
390     }
391
392     /* Now we have the sids that have to be shaken. We'll check the min and
393      * max atom numbers and this determines the shake block. DvdS 2007-07-19.
394      * For the purpose of making boundaries all atoms in between need to be
395      * part of the shake block too. There may be cases where blocks overlap
396      * and they will have to be merged.
397      */
398     merge_sid(0, at_end, nsid, sid, sblock);
399     sfree(sid);
400     /* Due to unknown reason this free generates a problem sometimes */
401     done_graph(g);
402     if (debug)
403     {
404         fprintf(debug, "Done gen_sblocks\n");
405     }
406 }