1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * GROningen Mixture of Alchemy and Childrens' Stories
54 typedef struct gmx_rmpbc {
62 static t_graph *gmx_rmpbc_get_graph(gmx_rmpbc_t gpbc, int ePBC, int natoms)
70 || gpbc->idef->ntypes <= 0)
76 for (i = 0; i < gpbc->ngraph; i++)
78 if (natoms == gpbc->graph[i].natoms)
85 /* We'd like to check with the number of atoms in the topology,
86 * but we don't have that available.
87 * So we check against the number of atoms that gmx_rmpbc_init
90 if (natoms > gpbc->natoms_init)
92 gmx_fatal(FARGS, "Structure or trajectory file has more atoms (%d) than the topology (%d)", natoms, gpbc->natoms_init);
95 srenew(gpbc->graph, gpbc->ngraph);
96 gr = &gpbc->graph[gpbc->ngraph-1];
98 gr->gr = mk_graph(NULL, gpbc->idef, 0, natoms, FALSE, FALSE);
104 gmx_rmpbc_t gmx_rmpbc_init(t_idef *idef, int ePBC, int natoms)
110 gpbc->natoms_init = natoms;
112 /* This sets pbc when we now it,
113 * otherwise we guess it from the instantaneous box in the trajectory.
118 if (gpbc->idef->ntypes <= 0)
122 "WARNING: If there are molecules in the input trajectory file\n"
123 " that are broken across periodic boundaries, they\n"
124 " cannot be made whole (or treated as whole) without\n"
125 " you providing a run input file.\n\n");
131 void gmx_rmpbc_done(gmx_rmpbc_t gpbc)
137 for (i = 0; i < gpbc->ngraph; i++)
139 done_graph(gpbc->graph[i].gr);
140 sfree(gpbc->graph[i].gr);
142 if (gpbc->graph != NULL)
150 static int gmx_rmpbc_ePBC(gmx_rmpbc_t gpbc, matrix box)
152 if (NULL != gpbc && gpbc->ePBC >= 0)
158 return guess_ePBC(box);
162 void gmx_rmpbc(gmx_rmpbc_t gpbc, int natoms, matrix box, rvec x[])
167 ePBC = gmx_rmpbc_ePBC(gpbc, box);
168 gr = gmx_rmpbc_get_graph(gpbc, ePBC, natoms);
171 mk_mshift(stdout, gr, ePBC, box, x);
172 shift_self(gr, box, x);
176 void gmx_rmpbc_copy(gmx_rmpbc_t gpbc, int natoms, matrix box, rvec x[], rvec x_s[])
182 ePBC = gmx_rmpbc_ePBC(gpbc, box);
183 gr = gmx_rmpbc_get_graph(gpbc, ePBC, natoms);
186 mk_mshift(stdout, gr, ePBC, box, x);
187 shift_x(gr, box, x, x_s);
191 for (i = 0; i < natoms; i++)
193 copy_rvec(x[i], x_s[i]);
198 void gmx_rmpbc_trxfr(gmx_rmpbc_t gpbc, t_trxframe *fr)
203 if (fr->bX && fr->bBox)
205 ePBC = gmx_rmpbc_ePBC(gpbc, fr->box);
206 gr = gmx_rmpbc_get_graph(gpbc, ePBC, fr->natoms);
209 mk_mshift(stdout, gr, ePBC, fr->box, fr->x);
210 shift_self(gr, fr->box, fr->x);
215 void rm_gropbc(t_atoms *atoms, rvec x[], matrix box)
220 /* check periodic boundary */
221 for (n = 1; (n < atoms->nr); n++)
223 for (m = DIM-1; m >= 0; m--)
225 dist = x[n][m]-x[n-1][m];
226 if (fabs(dist) > 0.9*box[m][m])
230 for (d = 0; d <= m; d++)
232 x[n][d] -= box[m][d];
237 for (d = 0; d <= m; d++)
239 x[n][d] += box[m][d];