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
,
111 gpbc
->natoms_init
= natoms
;
113 /* This sets pbc when we now it,
114 * otherwise we guess it from the instantaneous box in the trajectory.
119 if (gpbc
->idef
->ntypes
<= 0)
123 "WARNING: If there are molecules in the input trajectory file\n"
124 " that are broken across periodic boundaries, they\n"
125 " cannot be made whole (or treated as whole) without\n"
126 " you providing a run input file.\n\n");
132 void gmx_rmpbc_done(gmx_rmpbc_t gpbc
)
138 for(i
=0; i
<gpbc
->ngraph
; i
++)
140 done_graph(gpbc
->graph
[i
].gr
);
142 if (gpbc
->graph
!= NULL
)
149 static int gmx_rmpbc_ePBC(gmx_rmpbc_t gpbc
,matrix box
)
151 if (NULL
!= gpbc
&& gpbc
->ePBC
>= 0)
157 return guess_ePBC(box
);
161 void gmx_rmpbc(gmx_rmpbc_t gpbc
,int natoms
,matrix box
,rvec x
[])
166 ePBC
= gmx_rmpbc_ePBC(gpbc
,box
);
167 gr
= gmx_rmpbc_get_graph(gpbc
,ePBC
,natoms
);
170 mk_mshift(stdout
,gr
,ePBC
,box
,x
);
175 void gmx_rmpbc_copy(gmx_rmpbc_t gpbc
,int natoms
,matrix box
,rvec x
[],rvec x_s
[])
181 ePBC
= gmx_rmpbc_ePBC(gpbc
,box
);
182 gr
= gmx_rmpbc_get_graph(gpbc
,ePBC
,natoms
);
185 mk_mshift(stdout
,gr
,ePBC
,box
,x
);
186 shift_x(gr
,box
,x
,x_s
);
190 for(i
=0; i
<natoms
; i
++)
192 copy_rvec(x
[i
],x_s
[i
]);
197 void gmx_rmpbc_trxfr(gmx_rmpbc_t gpbc
,t_trxframe
*fr
)
202 if (fr
->bX
&& fr
->bBox
)
204 ePBC
= gmx_rmpbc_ePBC(gpbc
,fr
->box
);
205 gr
= gmx_rmpbc_get_graph(gpbc
,ePBC
,fr
->natoms
);
208 mk_mshift(stdout
,gr
,ePBC
,fr
->box
,fr
->x
);
209 shift_x(gr
,fr
->box
,fr
->x
,fr
->x
);
214 void rm_gropbc(t_atoms
*atoms
,rvec x
[],matrix box
)
219 /* check periodic boundary */
220 for(n
=1;(n
<atoms
->nr
);n
++)
222 for(m
=DIM
-1; m
>=0; m
--)
224 dist
= x
[n
][m
]-x
[n
-1][m
];
225 if (fabs(dist
) > 0.9*box
[m
][m
])
231 x
[n
][d
] -= box
[m
][d
];
238 x
[n
][d
] += box
[m
][d
];