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 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
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"
72 static t_graph
* gmx_rmpbc_get_graph(gmx_rmpbc_t gpbc
, int ePBC
, int natoms
)
77 if (ePBC
== epbcNONE
|| nullptr == gpbc
|| nullptr == gpbc
->idef
|| gpbc
->idef
->ntypes
<= 0)
83 for (i
= 0; i
< gpbc
->ngraph
; i
++)
85 if (natoms
== gpbc
->graph
[i
].natoms
)
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
97 if (natoms
> gpbc
->natoms_init
)
100 "Structure or trajectory file has more atoms (%d) than the topology (%d)",
101 natoms
, gpbc
->natoms_init
);
104 srenew(gpbc
->graph
, gpbc
->ngraph
);
105 gr
= &gpbc
->graph
[gpbc
->ngraph
- 1];
107 gr
->gr
= mk_graph(nullptr, gpbc
->idef
, 0, natoms
, FALSE
, FALSE
);
113 gmx_rmpbc_t
gmx_rmpbc_init(const t_idef
* idef
, int ePBC
, int natoms
)
119 gpbc
->natoms_init
= natoms
;
121 /* This sets pbc when we now it,
122 * otherwise we guess it from the instantaneous box in the trajectory.
127 if (gpbc
->idef
->ntypes
<= 0)
131 "WARNING: If there are molecules in the input trajectory file\n"
132 " that are broken across periodic boundaries, they\n"
133 " cannot be made whole (or treated as whole) without\n"
134 " you providing a run input file.\n\n");
140 void gmx_rmpbc_done(gmx_rmpbc_t gpbc
)
146 for (i
= 0; i
< gpbc
->ngraph
; i
++)
148 done_graph(gpbc
->graph
[i
].gr
);
149 sfree(gpbc
->graph
[i
].gr
);
151 if (gpbc
->graph
!= nullptr)
159 static int gmx_rmpbc_ePBC(gmx_rmpbc_t gpbc
, const matrix box
)
161 if (nullptr != gpbc
&& gpbc
->ePBC
>= 0)
167 return guess_ePBC(box
);
171 void gmx_rmpbc(gmx_rmpbc_t gpbc
, int natoms
, const matrix box
, rvec x
[])
176 ePBC
= gmx_rmpbc_ePBC(gpbc
, box
);
177 gr
= gmx_rmpbc_get_graph(gpbc
, ePBC
, natoms
);
180 mk_mshift(stdout
, gr
, ePBC
, box
, x
);
181 shift_self(gr
, box
, x
);
185 void gmx_rmpbc_copy(gmx_rmpbc_t gpbc
, int natoms
, const matrix box
, rvec x
[], rvec x_s
[])
191 ePBC
= gmx_rmpbc_ePBC(gpbc
, box
);
192 gr
= gmx_rmpbc_get_graph(gpbc
, ePBC
, natoms
);
195 mk_mshift(stdout
, gr
, ePBC
, box
, x
);
196 shift_x(gr
, box
, x
, x_s
);
200 for (i
= 0; i
< natoms
; i
++)
202 copy_rvec(x
[i
], x_s
[i
]);
207 void gmx_rmpbc_trxfr(gmx_rmpbc_t gpbc
, t_trxframe
* fr
)
212 if (fr
->bX
&& fr
->bBox
)
214 ePBC
= gmx_rmpbc_ePBC(gpbc
, fr
->box
);
215 gr
= gmx_rmpbc_get_graph(gpbc
, ePBC
, fr
->natoms
);
218 mk_mshift(stdout
, gr
, ePBC
, fr
->box
, fr
->x
);
219 shift_self(gr
, fr
->box
, fr
->x
);
224 void rm_gropbc(const t_atoms
* atoms
, rvec x
[], const matrix box
)
229 /* check periodic boundary */
230 for (n
= 1; (n
< atoms
->nr
); n
++)
232 for (m
= DIM
- 1; m
>= 0; m
--)
234 dist
= x
[n
][m
] - x
[n
- 1][m
];
235 if (std::abs(dist
) > 0.9 * box
[m
][m
])
239 for (d
= 0; d
<= m
; d
++)
241 x
[n
][d
] -= box
[m
][d
];
246 for (d
= 0; d
<= m
; d
++)
248 x
[n
][d
] += box
[m
][d
];