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