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, 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/pbcutil/mshift.h"
47 #include "gromacs/legacyheaders/types/ifunc.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/pbcutil/pbc.h"
51 #include "gromacs/utility/fatalerror.h"
52 #include "gromacs/utility/smalloc.h"
54 /************************************************************
56 * S H I F T U T I L I T I E S
58 ************************************************************/
61 /************************************************************
63 * G R A P H G E N E R A T I O N C O D E
65 ************************************************************/
67 static void add_gbond(t_graph *g, atom_id a0, atom_id a1)
73 inda0 = a0 - g->at_start;
74 inda1 = a1 - g->at_start;
76 /* Search for a direct edge between a0 and a1.
77 * All egdes are bidirectional, so we only need to search one way.
79 for (i = 0; (i < g->nedge[inda0] && !bFound); i++)
81 bFound = (g->edge[inda0][i] == a1);
86 g->edge[inda0][g->nedge[inda0]++] = a1;
87 g->edge[inda1][g->nedge[inda1]++] = a0;
91 static void mk_igraph(t_graph *g, int ftype, t_ilist *il,
92 int at_start, int at_end,
105 np = interaction_function[ftype].nratoms;
107 if (ia[1] >= at_start && ia[1] < at_end)
109 if (ia[np] >= at_end)
112 "Molecule in topology has atom numbers below and "
113 "above natoms (%d).\n"
114 "You are probably trying to use a trajectory which does "
115 "not match the first %d atoms of the run input file.\n"
116 "You can make a matching run input file with gmx convert-tpr.",
119 if (ftype == F_SETTLE)
121 /* Bond all the atoms in the settle */
122 add_gbond(g, ia[1], ia[2]);
123 add_gbond(g, ia[1], ia[3]);
125 else if (part == NULL)
127 /* Simply add this bond */
128 for (j = 1; j < np; j++)
130 add_gbond(g, ia[j], ia[j+1]);
135 /* Add this bond when it connects two unlinked parts of the graph */
136 for (j = 1; j < np; j++)
138 if (part[ia[j]] != part[ia[j+1]])
140 add_gbond(g, ia[j], ia[j+1]);
150 GMX_ATTRIBUTE_NORETURN static void g_error(int line, const char *file)
152 gmx_fatal(FARGS, "Tring to print non existant graph (file %s, line %d)",
155 #define GCHECK(g) if (g == NULL) g_error(__LINE__, __FILE__)
157 void p_graph(FILE *log, const char *title, t_graph *g)
160 const char *cc[egcolNR] = { "W", "G", "B" };
163 fprintf(log, "graph: %s\n", title);
164 fprintf(log, "nnodes: %d\n", g->nnodes);
165 fprintf(log, "nbound: %d\n", g->nbound);
166 fprintf(log, "start: %d\n", g->at_start);
167 fprintf(log, "end: %d\n", g->at_end);
168 fprintf(log, " atom shiftx shifty shiftz C nedg e1 e2 etc.\n");
169 for (i = 0; (i < g->nnodes); i++)
173 fprintf(log, "%5d%7d%7d%7d %1s%5d", g->at_start+i+1,
174 g->ishift[g->at_start+i][XX],
175 g->ishift[g->at_start+i][YY],
176 g->ishift[g->at_start+i][ZZ],
177 (g->negc > 0) ? cc[g->egc[i]] : " ",
179 for (j = 0; (j < g->nedge[i]); j++)
181 fprintf(log, " %5u", g->edge[i][j]+1);
189 static void calc_1se(t_graph *g, int ftype, t_ilist *il,
190 int nbond[], int at_start, int at_end)
192 int k, nratoms, end, j;
198 for (j = 0; (j < end); j += nratoms+1, ia += nratoms+1)
200 nratoms = interaction_function[ftype].nratoms;
202 if (ftype == F_SETTLE)
205 if (iaa >= at_start && iaa < at_end)
210 g->at_start = std::min(g->at_start, iaa);
211 g->at_end = std::max(g->at_end, iaa+2+1);
216 for (k = 1; (k <= nratoms); k++)
219 if (iaa >= at_start && iaa < at_end)
221 g->at_start = std::min(g->at_start, iaa);
222 g->at_end = std::max(g->at_end, iaa+1);
223 /* When making the graph we (might) link all atoms in an interaction
224 * sequentially. Therefore the end atoms add 1 to the count,
225 * the middle atoms 2.
227 if (k == 1 || k == nratoms)
241 static int calc_start_end(FILE *fplog, t_graph *g, t_ilist il[],
242 int at_start, int at_end,
247 g->at_start = at_end;
250 /* First add all the real bonds: they should determine the molecular
253 for (i = 0; (i < F_NRE); i++)
255 if (interaction_function[i].flags & IF_CHEMBOND)
257 calc_1se(g, i, &il[i], nbond, at_start, at_end);
260 /* Then add all the other interactions in fixed lists, but first
261 * check to see what's there already.
263 for (i = 0; (i < F_NRE); i++)
265 if (!(interaction_function[i].flags & IF_CHEMBOND))
267 calc_1se(g, i, &il[i], nbond, at_start, at_end);
273 for (i = g->at_start; (i < g->at_end); i++)
276 nnb = std::max(nnb, nbond[i]);
280 fprintf(fplog, "Max number of connections per atom is %d\n", nnb);
281 fprintf(fplog, "Total number of connections is %d\n", nbtot);
288 static void compact_graph(FILE *fplog, t_graph *g)
290 int i, j, n, max_nedge;
294 for (i = 0; i < g->nnodes; i++)
296 for (j = 0; j < g->nedge[i]; j++)
298 g->edge[0][n++] = g->edge[i][j];
300 max_nedge = std::max(max_nedge, g->nedge[i]);
302 srenew(g->edge[0], n);
303 /* set pointers after srenew because edge[0] might move */
304 for (i = 1; i < g->nnodes; i++)
306 g->edge[i] = g->edge[i-1] + g->nedge[i-1];
311 fprintf(fplog, "Max number of graph edges per atom is %d\n",
313 fprintf(fplog, "Total number of graph edges is %d\n", n);
317 static gmx_bool determine_graph_parts(t_graph *g, int *part)
321 atom_id at_i, *at_i2;
324 /* Initialize the part array with all entries different */
325 for (at_i = g->at_start; at_i < g->at_end; at_i++)
330 /* Loop over the graph until the part array is fixed */
335 for (i = 0; (i < g->nnodes); i++)
337 at_i = g->at_start + i;
339 for (e = 0; e < g->nedge[i]; e++)
341 /* Set part for both nodes to the minimum */
342 if (part[at_i2[e]] > part[at_i])
344 part[at_i2[e]] = part[at_i];
347 else if (part[at_i2[e]] < part[at_i])
349 part[at_i] = part[at_i2[e]];
353 if (part[at_i] != part[g->at_start])
360 fprintf(debug, "graph part[] nchanged=%d, bMultiPart=%d\n",
361 nchanged, bMultiPart);
364 while (nchanged > 0);
369 void mk_graph_ilist(FILE *fplog,
370 t_ilist *ilist, int at_start, int at_end,
371 gmx_bool bShakeOnly, gmx_bool bSettle,
378 /* The naming is somewhat confusing, but we need g->at0 and g->at1
379 * for shifthing coordinates to a new array (not in place) when
380 * some atoms are not connected by the graph, which runs from
381 * g->at_start (>= g->at0) to g->at_end (<= g->at1).
387 nbtot = calc_start_end(fplog, g, ilist, at_start, at_end, nbond);
389 if (g->at_start >= g->at_end)
391 g->at_start = at_start;
398 g->nnodes = g->at_end - g->at_start;
399 snew(g->nedge, g->nnodes);
400 snew(g->edge, g->nnodes);
401 /* Allocate a single array and set pointers into it */
402 snew(g->edge[0], nbtot);
403 for (i = 1; (i < g->nnodes); i++)
405 g->edge[i] = g->edge[i-1] + nbond[g->at_start+i-1];
410 /* First add all the real bonds: they should determine the molecular
413 for (i = 0; (i < F_NRE); i++)
415 if (interaction_function[i].flags & IF_CHEMBOND)
417 mk_igraph(g, i, &(ilist[i]), at_start, at_end, NULL);
421 /* Determine of which separated parts the IF_CHEMBOND graph consists.
422 * Store the parts in the nbond array.
424 bMultiPart = determine_graph_parts(g, nbond);
428 /* Then add all the other interactions in fixed lists,
429 * but only when they connect parts of the graph
430 * that are not connected through IF_CHEMBOND interactions.
432 for (i = 0; (i < F_NRE); i++)
434 if (!(interaction_function[i].flags & IF_CHEMBOND))
436 mk_igraph(g, i, &(ilist[i]), at_start, at_end, nbond);
441 /* Removed all the unused space from the edge array */
442 compact_graph(fplog, g);
446 /* This is a special thing used in splitter.c to generate shake-blocks */
447 mk_igraph(g, F_CONSTR, &(ilist[F_CONSTR]), at_start, at_end, NULL);
450 mk_igraph(g, F_SETTLE, &(ilist[F_SETTLE]), at_start, at_end, NULL);
454 for (i = 0; (i < g->nnodes); i++)
468 snew(g->ishift, g->at1);
472 p_graph(debug, "graph", g);
476 t_graph *mk_graph(FILE *fplog,
477 t_idef *idef, int at_start, int at_end,
478 gmx_bool bShakeOnly, gmx_bool bSettle)
484 mk_graph_ilist(fplog, idef->il, at_start, at_end, bShakeOnly, bSettle, g);
489 void done_graph(t_graph *g)
502 /************************************************************
504 * S H I F T C A L C U L A T I O N C O D E
506 ************************************************************/
508 static void mk_1shift_tric(int npbcdim, matrix box, rvec hbox,
509 rvec xi, rvec xj, int *mi, int *mj)
511 /* Calculate periodicity for triclinic box... */
515 rvec_sub(xi, xj, dx);
518 for (m = npbcdim-1; (m >= 0); m--)
520 /* If dx < hbox, then xj will be reduced by box, so that
521 * xi - xj will be bigger
523 if (dx[m] < -hbox[m])
526 for (d = m-1; d >= 0; d--)
531 else if (dx[m] >= hbox[m])
534 for (d = m-1; d >= 0; d--)
546 static void mk_1shift(int npbcdim, rvec hbox, rvec xi, rvec xj, int *mi, int *mj)
548 /* Calculate periodicity for rectangular box... */
552 rvec_sub(xi, xj, dx);
555 for (m = 0; (m < npbcdim); m++)
557 /* If dx < hbox, then xj will be reduced by box, so that
558 * xi - xj will be bigger
560 if (dx[m] < -hbox[m])
564 else if (dx[m] >= hbox[m])
575 static void mk_1shift_screw(matrix box, rvec hbox,
576 rvec xi, rvec xj, int *mi, int *mj)
578 /* Calculate periodicity for rectangular box... */
582 if ((mi[XX] > 0 && mi[XX] % 2 == 1) ||
583 (mi[XX] < 0 && -mi[XX] % 2 == 1))
592 rvec_sub(xi, xj, dx);
594 if (dx[XX] < -hbox[XX])
598 else if (dx[XX] >= hbox[XX])
606 if (mj[XX] != mi[XX])
609 dx[YY] = xi[YY] - (box[YY][YY] + box[ZZ][YY] - xj[YY]);
610 dx[ZZ] = xi[ZZ] - (box[ZZ][ZZ] - xj[ZZ]);
612 for (m = 1; (m < DIM); m++)
614 /* The signs are taken such that we can first shift x and rotate
615 * and then shift y and z.
617 if (dx[m] < -hbox[m])
619 mj[m] = mi[m] - signi;
621 else if (dx[m] >= hbox[m])
623 mj[m] = mi[m] + signi;
632 static int mk_grey(egCol egc[], t_graph *g, int *AtomI,
633 int npbcdim, matrix box, rvec x[], int *nerror)
635 int m, j, ng, ai, aj, g0;
641 for (m = 0; (m < DIM); m++)
643 hbox[m] = box[m][m]*0.5;
645 bTriclinic = TRICLINIC(box);
651 /* Loop over all the bonds */
652 for (j = 0; (j < g->nedge[ai-g0]); j++)
654 aj = g->edge[ai-g0][j];
655 /* If there is a white one, make it grey and set pbc */
658 mk_1shift_screw(box, hbox, x[ai], x[aj], g->ishift[ai], is_aj);
662 mk_1shift_tric(npbcdim, box, hbox, x[ai], x[aj], g->ishift[ai], is_aj);
666 mk_1shift(npbcdim, hbox, x[ai], x[aj], g->ishift[ai], is_aj);
669 if (egc[aj-g0] == egcolWhite)
671 if (aj - g0 < *AtomI)
675 egc[aj-g0] = egcolGrey;
677 copy_ivec(is_aj, g->ishift[aj]);
681 else if ((is_aj[XX] != g->ishift[aj][XX]) ||
682 (is_aj[YY] != g->ishift[aj][YY]) ||
683 (is_aj[ZZ] != g->ishift[aj][ZZ]))
687 set_pbc(&pbc, -1, box);
688 pbc_dx(&pbc, x[ai], x[aj], dx);
689 fprintf(debug, "mk_grey: shifts for atom %d due to atom %d\n"
690 "are (%d,%d,%d), should be (%d,%d,%d)\n"
692 aj+1, ai+1, is_aj[XX], is_aj[YY], is_aj[ZZ],
693 g->ishift[aj][XX], g->ishift[aj][YY], g->ishift[aj][ZZ],
694 dx[XX], dx[YY], dx[ZZ]);
702 static int first_colour(int fC, egCol Col, t_graph *g, egCol egc[])
703 /* Return the first node with colour Col starting at fC.
704 * return -1 if none found.
709 for (i = fC; (i < g->nnodes); i++)
711 if ((g->nedge[i] > 0) && (egc[i] == Col))
720 void mk_mshift(FILE *log, t_graph *g, int ePBC, matrix box, rvec x[])
722 static int nerror_tot = 0;
725 int nW, nG, nB; /* Number of Grey, Black, White */
726 int fW, fG; /* First of each category */
729 g->bScrewPBC = (ePBC == epbcSCREW);
741 /* This puts everything in the central box, that is does not move it
742 * at all. If we return without doing this for a system without bonds
743 * (i.e. only settles) all water molecules are moved to the opposite octant
745 for (i = g->at0; (i < g->at1); i++)
747 g->ishift[i][XX] = g->ishift[i][YY] = g->ishift[i][ZZ] = 0;
756 if (nnodes > g->negc)
759 srenew(g->egc, g->negc);
761 memset(g->egc, 0, (size_t)(nnodes*sizeof(g->egc[0])));
769 /* We even have a loop invariant:
770 * nW+nG+nB == g->nbound
773 fprintf(log, "Starting W loop\n");
777 /* Find the first white, this will allways be a larger
778 * number than before, because no nodes are made white
781 if ((fW = first_colour(fW, egcolWhite, g, g->egc)) == -1)
783 gmx_fatal(FARGS, "No WHITE nodes found while nW=%d\n", nW);
786 /* Make the first white node grey */
787 g->egc[fW] = egcolGrey;
791 /* Initial value for the first grey */
794 fprintf(log, "Starting G loop (nW=%d, nG=%d, nB=%d, total %d)\n",
795 nW, nG, nB, nW+nG+nB);
799 if ((fG = first_colour(fG, egcolGrey, g, g->egc)) == -1)
801 gmx_fatal(FARGS, "No GREY nodes found while nG=%d\n", nG);
804 /* Make the first grey node black */
805 g->egc[fG] = egcolBlack;
809 /* Make all the neighbours of this black node grey
810 * and set their periodicity
812 ng = mk_grey(g->egc, g, &fG, npbcdim, box, x, &nerror);
813 /* ng is the number of white nodes made grey */
821 if (nerror_tot <= 100)
823 fprintf(stderr, "There were %d inconsistent shifts. Check your topology\n",
827 fprintf(log, "There were %d inconsistent shifts. Check your topology\n",
831 if (nerror_tot == 100)
833 fprintf(stderr, "Will stop reporting inconsistent shifts\n");
836 fprintf(log, "Will stop reporting inconsistent shifts\n");
842 /************************************************************
844 * A C T U A L S H I F T C O D E
846 ************************************************************/
848 void shift_x(t_graph *g, matrix box, rvec x[], rvec x_s[])
859 for (j = g->at0; j < g0; j++)
861 copy_rvec(x[j], x_s[j]);
866 for (j = g0; (j < g1); j++)
872 if ((tx > 0 && tx % 2 == 1) ||
873 (tx < 0 && -tx %2 == 1))
875 x_s[j][XX] = x[j][XX] + tx*box[XX][XX];
876 x_s[j][YY] = box[YY][YY] + box[ZZ][YY] - x[j][YY];
877 x_s[j][ZZ] = box[ZZ][ZZ] - x[j][ZZ];
881 x_s[j][XX] = x[j][XX];
883 x_s[j][YY] = x[j][YY] + ty*box[YY][YY] + tz*box[ZZ][YY];
884 x_s[j][ZZ] = x[j][ZZ] + tz*box[ZZ][ZZ];
887 else if (TRICLINIC(box))
889 for (j = g0; (j < g1); j++)
895 x_s[j][XX] = x[j][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
896 x_s[j][YY] = x[j][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
897 x_s[j][ZZ] = x[j][ZZ]+tz*box[ZZ][ZZ];
902 for (j = g0; (j < g1); j++)
908 x_s[j][XX] = x[j][XX]+tx*box[XX][XX];
909 x_s[j][YY] = x[j][YY]+ty*box[YY][YY];
910 x_s[j][ZZ] = x[j][ZZ]+tz*box[ZZ][ZZ];
914 for (j = g1; j < g->at1; j++)
916 copy_rvec(x[j], x_s[j]);
920 void shift_self(t_graph *g, matrix box, rvec x[])
928 gmx_incons("screw pbc not implemented for shift_self");
936 fprintf(stderr, "Shifting atoms %d to %d\n", g0, g0+gn);
940 for (j = g0; (j < g1); j++)
946 x[j][XX] = x[j][XX]+tx*box[XX][XX]+ty*box[YY][XX]+tz*box[ZZ][XX];
947 x[j][YY] = x[j][YY]+ty*box[YY][YY]+tz*box[ZZ][YY];
948 x[j][ZZ] = x[j][ZZ]+tz*box[ZZ][ZZ];
953 for (j = g0; (j < g1); j++)
959 x[j][XX] = x[j][XX]+tx*box[XX][XX];
960 x[j][YY] = x[j][YY]+ty*box[YY][YY];
961 x[j][ZZ] = x[j][ZZ]+tz*box[ZZ][ZZ];
966 void unshift_x(t_graph *g, matrix box, rvec x[], rvec x_s[])
974 gmx_incons("screw pbc not implemented (yet) for unshift_x");
981 for (j = g->at0; j < g0; j++)
983 copy_rvec(x_s[j], x[j]);
988 for (j = g0; (j < g1); j++)
994 x[j][XX] = x_s[j][XX]-tx*box[XX][XX]-ty*box[YY][XX]-tz*box[ZZ][XX];
995 x[j][YY] = x_s[j][YY]-ty*box[YY][YY]-tz*box[ZZ][YY];
996 x[j][ZZ] = x_s[j][ZZ]-tz*box[ZZ][ZZ];
1001 for (j = g0; (j < g1); j++)
1007 x[j][XX] = x_s[j][XX]-tx*box[XX][XX];
1008 x[j][YY] = x_s[j][YY]-ty*box[YY][YY];
1009 x[j][ZZ] = x_s[j][ZZ]-tz*box[ZZ][ZZ];
1013 for (j = g1; j < g->at1; j++)
1015 copy_rvec(x_s[j], x[j]);
1019 void unshift_self(t_graph *g, matrix box, rvec x[])
1027 gmx_incons("screw pbc not implemented for unshift_self");
1036 for (j = g0; (j < g1); j++)
1042 x[j][XX] = x[j][XX]-tx*box[XX][XX]-ty*box[YY][XX]-tz*box[ZZ][XX];
1043 x[j][YY] = x[j][YY]-ty*box[YY][YY]-tz*box[ZZ][YY];
1044 x[j][ZZ] = x[j][ZZ]-tz*box[ZZ][ZZ];
1049 for (j = g0; (j < g1); j++)
1055 x[j][XX] = x[j][XX]-tx*box[XX][XX];
1056 x[j][YY] = x[j][YY]-ty*box[YY][YY];
1057 x[j][ZZ] = x[j][ZZ]-tz*box[ZZ][ZZ];