2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
39 #include "gromacs/legacyheaders/splitter.h"
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"
58 static int sid_comp(const void *a, const void *b)
69 return (sa->atom-sb->atom);
77 static int mk_grey(egCol egc[], t_graph *g, int *AtomI,
78 int maxsid, t_sid sid[])
80 int j, ng, ai, aj, g0;
86 /* Loop over all the bonds */
87 for (j = 0; (j < g->nedge[ai]); j++)
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)
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)
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__);
109 sid[aj+g0].sid = sid[ai+g0].sid;
110 sid[aj+g0].atom = aj+g0;
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.
125 for (i = fC; (i < g->nnodes); i++)
127 if ((g->nedge[i] > 0) && (egc[i] == Col))
136 static int mk_sblocks(FILE *fp, t_graph *g, int maxsid, t_sid sid[])
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 */
161 /* We even have a loop invariant:
162 * nW+nG+nB == g->nbound
167 fprintf(fp, "Walking down the molecule graph to make constraint-blocks\n");
172 /* Find the first white, this will allways be a larger
173 * number than before, because no nodes are made white
176 if ((fW = first_colour(fW, egcolWhite, g, egc)) == -1)
178 gmx_fatal(FARGS, "No WHITE nodes found while nW=%d\n", nW);
181 /* Make the first white node grey, and set the block number */
183 range_check(fW+g0, 0, maxsid);
184 sid[fW+g0].sid = nblock++;
188 /* Initial value for the first grey */
193 fprintf(debug, "Starting G loop (nW=%d, nG=%d, nB=%d, total %d)\n",
194 nW, nG, nB, nW+nG+nB);
199 if ((fG = first_colour(fG, egcolGrey, g, egc)) == -1)
201 gmx_fatal(FARGS, "No GREY nodes found while nG=%d\n", nG);
204 /* Make the first grey node black */
205 egc[fG] = egcolBlack;
209 /* Make all the neighbours of this black node grey
210 * and set their block number
212 ng = mk_grey(egc, g, &fG, maxsid, sid);
213 /* ng is the number of white nodes made grey */
222 fprintf(debug, "Found %d shake blocks\n", nblock);
230 int first, last, sid;
233 static int ms_comp(const void *a, const void *b)
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);
239 d = ma->first-mb->first;
242 return ma->last-mb->last;
250 static int merge_sid(int at_start, int at_end, int nsid, t_sid sid[],
253 int i, j, k, n, isid, ndel;
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
261 /* Determine first and last atom in each shake ID */
264 for (k = 0; (k < nsid); k++)
266 ms[k].first = at_end+1;
270 for (i = at_start; (i < at_end); i++)
273 range_check(isid, -1, nsid);
276 ms[isid].first = std::min(ms[isid].first, sid[i].atom);
277 ms[isid].last = std::max(ms[isid].last, sid[i].atom);
280 qsort(ms, nsid, sizeof(ms[0]), ms_comp);
282 /* Now merge the overlapping ones */
284 for (k = 0; (k < nsid); )
286 for (j = k+1; (j < nsid); )
288 if (ms[j].first <= ms[k].last)
290 ms[k].last = std::max(ms[k].last, ms[j].last);
291 ms[k].first = std::min(ms[k].first, ms[j].first);
307 for (k = 0; (k < nsid); k++)
309 while ((k < nsid-1) && (ms[k].sid == -1))
311 for (j = k+1; (j < nsid); j++)
313 std::memcpy(&(ms[j-1]), &(ms[j]), sizeof(ms[0]));
319 for (k = at_start; (k < at_end); k++)
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++)
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++)
336 range_check(n, 0, sblock->nra);
338 range_check(j, 0, at_end);
339 if (sid[j].sid == -1)
345 fprintf(stderr, "Double sids (%d, %d) for atom %d\n", sid[j].sid, k, j);
351 GMX_RELEASE_ASSERT(sblock->index[k] == sblock->nra, "Internal inconsistency; sid merge failed");
358 void gen_sblocks(FILE *fp, int at_start, int at_end,
359 t_idef *idef, t_blocka *sblock,
367 g = mk_graph(NULL, idef, at_start, at_end, TRUE, bSettle);
370 p_graph(debug, "Graaf Dracula", g);
373 for (i = at_start; (i < at_end); i++)
378 nsid = mk_sblocks(fp, g, at_end, sid);
385 /* Now sort the shake blocks... */
386 qsort(sid+at_start, at_end-at_start, static_cast<size_t>(sizeof(sid[0])), sid_comp);
390 fprintf(debug, "Sorted shake block\n");
391 for (i = at_start; (i < at_end); i++)
393 fprintf(debug, "sid[%5d] = atom:%5d sid:%5d\n", i, sid[i].atom, sid[i].sid);
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++)
399 if (sid[i0].sid > -1)
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.
411 merge_sid(at_start, at_end, nsid, sid, sblock);
413 /* Due to unknown reason this free generates a problem sometimes */
418 fprintf(debug, "Done gen_sblocks\n");