Create separate module for trajectory data
[gromacs.git] / src / gromacs / pbcutil / rmpbc.cpp
blob4440ec87a41f36b7760cd2912804aa0f7cebf241
1 /*
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.
37 #include "gmxpre.h"
39 #include "rmpbc.h"
41 #include <stdlib.h>
43 #include <cmath>
45 #include <algorithm>
47 #include "gromacs/math/vec.h"
48 #include "gromacs/pbcutil/mshift.h"
49 #include "gromacs/pbcutil/pbc.h"
50 #include "gromacs/topology/atoms.h"
51 #include "gromacs/topology/idef.h"
52 #include "gromacs/trajectory/trajectoryframe.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/futil.h"
55 #include "gromacs/utility/smalloc.h"
57 typedef struct {
58 int natoms;
59 t_graph *gr;
60 } rmpbc_graph_t;
62 struct gmx_rmpbc {
63 t_idef *idef;
64 int natoms_init;
65 int ePBC;
66 int ngraph;
67 rmpbc_graph_t *graph;
70 static t_graph *gmx_rmpbc_get_graph(gmx_rmpbc_t gpbc, int ePBC, int natoms)
72 int i;
73 rmpbc_graph_t *gr;
75 if (ePBC == epbcNONE
76 || NULL == gpbc
77 || NULL == gpbc->idef
78 || gpbc->idef->ntypes <= 0)
80 return NULL;
83 gr = NULL;
84 for (i = 0; i < gpbc->ngraph; i++)
86 if (natoms == gpbc->graph[i].natoms)
88 gr = &gpbc->graph[i];
91 if (gr == NULL)
93 /* We'd like to check with the number of atoms in the topology,
94 * but we don't have that available.
95 * So we check against the number of atoms that gmx_rmpbc_init
96 * was called with.
98 if (natoms > gpbc->natoms_init)
100 gmx_fatal(FARGS, "Structure or trajectory file has more atoms (%d) than the topology (%d)", natoms, gpbc->natoms_init);
102 gpbc->ngraph++;
103 srenew(gpbc->graph, gpbc->ngraph);
104 gr = &gpbc->graph[gpbc->ngraph-1];
105 gr->natoms = natoms;
106 gr->gr = mk_graph(NULL, gpbc->idef, 0, natoms, FALSE, FALSE);
109 return gr->gr;
112 gmx_rmpbc_t gmx_rmpbc_init(t_idef *idef, int ePBC, int natoms)
114 gmx_rmpbc_t gpbc;
116 snew(gpbc, 1);
118 gpbc->natoms_init = natoms;
120 /* This sets pbc when we now it,
121 * otherwise we guess it from the instantaneous box in the trajectory.
123 gpbc->ePBC = ePBC;
125 gpbc->idef = idef;
126 if (gpbc->idef->ntypes <= 0)
128 fprintf(stderr,
129 "\n"
130 "WARNING: If there are molecules in the input trajectory file\n"
131 " that are broken across periodic boundaries, they\n"
132 " cannot be made whole (or treated as whole) without\n"
133 " you providing a run input file.\n\n");
136 return gpbc;
139 void gmx_rmpbc_done(gmx_rmpbc_t gpbc)
141 int i;
143 if (NULL != gpbc)
145 for (i = 0; i < gpbc->ngraph; i++)
147 done_graph(gpbc->graph[i].gr);
148 sfree(gpbc->graph[i].gr);
150 if (gpbc->graph != NULL)
152 sfree(gpbc->graph);
154 sfree(gpbc);
158 static int gmx_rmpbc_ePBC(gmx_rmpbc_t gpbc, matrix box)
160 if (NULL != gpbc && gpbc->ePBC >= 0)
162 return gpbc->ePBC;
164 else
166 return guess_ePBC(box);
170 void gmx_rmpbc(gmx_rmpbc_t gpbc, int natoms, matrix box, rvec x[])
172 int ePBC;
173 t_graph *gr;
175 ePBC = gmx_rmpbc_ePBC(gpbc, box);
176 gr = gmx_rmpbc_get_graph(gpbc, ePBC, natoms);
177 if (gr != NULL)
179 mk_mshift(stdout, gr, ePBC, box, x);
180 shift_self(gr, box, x);
184 void gmx_rmpbc_copy(gmx_rmpbc_t gpbc, int natoms, matrix box, rvec x[], rvec x_s[])
186 int ePBC;
187 t_graph *gr;
188 int i;
190 ePBC = gmx_rmpbc_ePBC(gpbc, box);
191 gr = gmx_rmpbc_get_graph(gpbc, ePBC, natoms);
192 if (gr != NULL)
194 mk_mshift(stdout, gr, ePBC, box, x);
195 shift_x(gr, box, x, x_s);
197 else
199 for (i = 0; i < natoms; i++)
201 copy_rvec(x[i], x_s[i]);
206 void gmx_rmpbc_trxfr(gmx_rmpbc_t gpbc, t_trxframe *fr)
208 int ePBC;
209 t_graph *gr;
211 if (fr->bX && fr->bBox)
213 ePBC = gmx_rmpbc_ePBC(gpbc, fr->box);
214 gr = gmx_rmpbc_get_graph(gpbc, ePBC, fr->natoms);
215 if (gr != NULL)
217 mk_mshift(stdout, gr, ePBC, fr->box, fr->x);
218 shift_self(gr, fr->box, fr->x);
223 void rm_gropbc(t_atoms *atoms, rvec x[], matrix box)
225 real dist;
226 int n, m, d;
228 /* check periodic boundary */
229 for (n = 1; (n < atoms->nr); n++)
231 for (m = DIM-1; m >= 0; m--)
233 dist = x[n][m]-x[n-1][m];
234 if (std::abs(dist) > 0.9*box[m][m])
236 if (dist > 0)
238 for (d = 0; d <= m; d++)
240 x[n][d] -= box[m][d];
243 else
245 for (d = 0; d <= m; d++)
247 x[n][d] += box[m][d];