Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / pbcutil / rmpbc.cpp
blob20336c6c73508ae7beec88a62a76064141d8aafe
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,2016,2017,2018, 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 <cmath>
42 #include <cstdlib>
44 #include <algorithm>
46 #include "gromacs/math/vec.h"
47 #include "gromacs/pbcutil/mshift.h"
48 #include "gromacs/pbcutil/pbc.h"
49 #include "gromacs/topology/atoms.h"
50 #include "gromacs/topology/idef.h"
51 #include "gromacs/trajectory/trajectoryframe.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/futil.h"
54 #include "gromacs/utility/smalloc.h"
56 typedef struct {
57 int natoms;
58 t_graph *gr;
59 } rmpbc_graph_t;
61 struct gmx_rmpbc {
62 const t_idef *idef;
63 int natoms_init;
64 int ePBC;
65 int ngraph;
66 rmpbc_graph_t *graph;
69 static t_graph *gmx_rmpbc_get_graph(gmx_rmpbc_t gpbc, int ePBC, int natoms)
71 int i;
72 rmpbc_graph_t *gr;
74 if (ePBC == epbcNONE
75 || nullptr == gpbc
76 || nullptr == gpbc->idef
77 || gpbc->idef->ntypes <= 0)
79 return nullptr;
82 gr = nullptr;
83 for (i = 0; i < gpbc->ngraph; i++)
85 if (natoms == gpbc->graph[i].natoms)
87 gr = &gpbc->graph[i];
90 if (gr == nullptr)
92 /* We'd like to check with the number of atoms in the topology,
93 * but we don't have that available.
94 * So we check against the number of atoms that gmx_rmpbc_init
95 * was called with.
97 if (natoms > gpbc->natoms_init)
99 gmx_fatal(FARGS, "Structure or trajectory file has more atoms (%d) than the topology (%d)", natoms, gpbc->natoms_init);
101 gpbc->ngraph++;
102 srenew(gpbc->graph, gpbc->ngraph);
103 gr = &gpbc->graph[gpbc->ngraph-1];
104 gr->natoms = natoms;
105 gr->gr = mk_graph(nullptr, gpbc->idef, 0, natoms, FALSE, FALSE);
108 return gr->gr;
111 gmx_rmpbc_t gmx_rmpbc_init(const t_idef *idef, int ePBC, int natoms)
113 gmx_rmpbc_t gpbc;
115 snew(gpbc, 1);
117 gpbc->natoms_init = natoms;
119 /* This sets pbc when we now it,
120 * otherwise we guess it from the instantaneous box in the trajectory.
122 gpbc->ePBC = ePBC;
124 gpbc->idef = idef;
125 if (gpbc->idef->ntypes <= 0)
127 fprintf(stderr,
128 "\n"
129 "WARNING: If there are molecules in the input trajectory file\n"
130 " that are broken across periodic boundaries, they\n"
131 " cannot be made whole (or treated as whole) without\n"
132 " you providing a run input file.\n\n");
135 return gpbc;
138 void gmx_rmpbc_done(gmx_rmpbc_t gpbc)
140 int i;
142 if (nullptr != gpbc)
144 for (i = 0; i < gpbc->ngraph; i++)
146 done_graph(gpbc->graph[i].gr);
147 sfree(gpbc->graph[i].gr);
149 if (gpbc->graph != nullptr)
151 sfree(gpbc->graph);
153 sfree(gpbc);
157 static int gmx_rmpbc_ePBC(gmx_rmpbc_t gpbc, const matrix box)
159 if (nullptr != gpbc && gpbc->ePBC >= 0)
161 return gpbc->ePBC;
163 else
165 return guess_ePBC(box);
169 void gmx_rmpbc(gmx_rmpbc_t gpbc, int natoms, const matrix box, rvec x[])
171 int ePBC;
172 t_graph *gr;
174 ePBC = gmx_rmpbc_ePBC(gpbc, box);
175 gr = gmx_rmpbc_get_graph(gpbc, ePBC, natoms);
176 if (gr != nullptr)
178 mk_mshift(stdout, gr, ePBC, box, x);
179 shift_self(gr, box, x);
183 void gmx_rmpbc_copy(gmx_rmpbc_t gpbc, int natoms, const matrix box, rvec x[], rvec x_s[])
185 int ePBC;
186 t_graph *gr;
187 int i;
189 ePBC = gmx_rmpbc_ePBC(gpbc, box);
190 gr = gmx_rmpbc_get_graph(gpbc, ePBC, natoms);
191 if (gr != nullptr)
193 mk_mshift(stdout, gr, ePBC, box, x);
194 shift_x(gr, box, x, x_s);
196 else
198 for (i = 0; i < natoms; i++)
200 copy_rvec(x[i], x_s[i]);
205 void gmx_rmpbc_trxfr(gmx_rmpbc_t gpbc, t_trxframe *fr)
207 int ePBC;
208 t_graph *gr;
210 if (fr->bX && fr->bBox)
212 ePBC = gmx_rmpbc_ePBC(gpbc, fr->box);
213 gr = gmx_rmpbc_get_graph(gpbc, ePBC, fr->natoms);
214 if (gr != nullptr)
216 mk_mshift(stdout, gr, ePBC, fr->box, fr->x);
217 shift_self(gr, fr->box, fr->x);
222 void rm_gropbc(const t_atoms *atoms, rvec x[], const matrix box)
224 real dist;
225 int n, m, d;
227 /* check periodic boundary */
228 for (n = 1; (n < atoms->nr); n++)
230 for (m = DIM-1; m >= 0; m--)
232 dist = x[n][m]-x[n-1][m];
233 if (std::abs(dist) > 0.9*box[m][m])
235 if (dist > 0)
237 for (d = 0; d <= m; d++)
239 x[n][d] -= box[m][d];
242 else
244 for (d = 0; d <= m; d++)
246 x[n][d] += box[m][d];