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) 2012,2014,2015,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020,2021, 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.
38 #ifndef GMX_PBCUTIL_PBC_H
39 #define GMX_PBCUTIL_PBC_H
45 #include "gromacs/math/vectypes.h"
46 #include "gromacs/utility/basedefinitions.h"
47 #include "gromacs/utility/enumerationhelpers.h"
48 #include "gromacs/utility/real.h"
59 //! Enumeration that contains all supported periodic boundary setups.
60 enum class PbcType : int
62 Xyz = 0, //!< Periodic boundaries in all dimensions.
63 No = 1, //!< No periodic boundaries.
64 XY = 2, //!< Only two dimensions are periodic.
65 Screw = 3, //!< Screw.
66 Unset = 4, //!< The type of PBC is not set or invalid.
70 //! Names for all values in PBC types enumeration
71 extern const gmx::EnumerationArray<PbcType, std::string> c_pbcTypeNames;
73 /* Maximum number of combinations of single triclinic box vectors
74 * required to shift atoms that are within a brick of the size of
75 * the diagonal of the box to within the maximum cut-off distance.
77 #define MAX_NTRICVEC 12
79 /*! \brief Structure containing info on periodic boundary conditions */
84 //! Number of dimensions in which PBC is exerted
86 /*! \brief Determines how to compute distance vectors.
88 * Indicator of how to compute distance vectors, depending
89 * on PBC type (depends on pbcType and dimensions with(out) DD)
93 /*! \brief Used for selecting which dimensions to use in PBC.
95 * In case of 1-D PBC this indicates which dimension is used,
96 * in case of 2-D PBC this indicates the opposite
99 //! The simulation box
101 //! The lengths of the diagonal of the full box
103 //! Halve of the above
105 //! Negative of the above
107 //! Maximum allowed cutoff squared for the box and PBC used
109 /*! \brief Number of triclinic shift vectors.
111 * Number of triclinic shift vectors depends on the skewedness
112 * of the box, that is mostly on the angles. For triclinic boxes
113 * we first take the closest image along each Cartesian dimension
114 * independently. When the resulting distance^2 is larger than
115 * max_cutoff2, up to ntric_vec triclinic shift vectors need to
116 * be tried. Because of the restrictions imposed on the unit-cell
117 * by GROMACS, ntric_vec <= MAX_NTRICVEC = 12.
120 //! The triclinic shift vectors in grid cells. Internal use only.
121 ivec tric_shift[MAX_NTRICVEC];
122 //! The triclinic shift vectors in length units
123 rvec tric_vec[MAX_NTRICVEC];
126 #define TRICLINIC(box) ((box)[YY][XX] != 0 || (box)[ZZ][XX] != 0 || (box)[ZZ][YY] != 0)
134 ecenterTRIC, /* 0.5*(a+b+c) */
135 ecenterRECT, /* (0.5*a[x],0.5*b[y],0.5*c[z]) */
136 ecenterZERO, /* (0,0,0) */
137 ecenterDEF = ecenterTRIC
140 /*! \brief Returns the number of dimensions that use pbc
142 * \param[in] pbcType The periodic boundary condition type
143 * \return the number of dimensions that use pbc, starting at X
145 int numPbcDimensions(PbcType pbcType);
147 /*! \brief Dump the contents of the pbc structure to the file
149 * \param[in] fp The file pointer to write to
150 * \param[in] pbc The periodic boundary condition information structure
152 void dump_pbc(FILE* fp, t_pbc* pbc);
154 /*! \brief Check the box for consistency
156 * When \p pbcType=PbcTypes::Unset, the type of pbc is guessed from the box matrix.
158 * \param[in] pbcType The pbc identifier
159 * \param[in] box The box matrix
160 * \return NULL if the box is supported by Gromacs.
161 * Otherwise returns a string with the problem.
163 const char* check_box(PbcType pbcType, const matrix box);
165 /*! \brief Creates box matrix from edge lengths and angles.
167 * \param[in,out] box The box matrix
168 * \param[in] vec The edge lengths
169 * \param[in] angleInDegrees The angles
171 void matrix_convert(matrix box, const rvec vec, const rvec angleInDegrees);
173 /*! \brief Compute the maximum cutoff for the box
175 * Returns the square of the maximum cut-off allowed for the box,
176 * taking into account that the grid neighborsearch code and pbc_dx
177 * only check combinations of single box-vector shifts.
179 * \param[in] pbcType The pbc identifier
180 * \param[in] box The box matrix
181 * \return the maximum cut-off.
183 real max_cutoff2(PbcType pbcType, const matrix box);
185 /*! \brief Guess PBC type
187 * Guesses the type of periodic boundary conditions using the box
189 * \param[in] box The box matrix
190 * \return The pbc type identifier
192 PbcType guessPbcType(const matrix box);
194 /*! \brief Corrects the box if necessary
196 * Checks for un-allowed box angles and corrects the box.
198 * \param[in] fplog File for debug output
199 * \param[in] step The MD step number
200 * \param[in] box The simulation cell
201 * \return TRUE when the box was corrected.
203 bool correct_box(FILE* fplog, int step, tensor box);
205 /*! \brief Initiate the periodic boundary condition algorithms.
207 * pbc_dx will not use pbc and return the normal difference vector
208 * when one or more of the diagonal elements of box are zero.
209 * When \p pbcType=PbcType::Unset, the type of pbc is guessed from the box matrix.
211 * \param[in,out] pbc The pbc information structure
212 * \param[in] pbcType The PBC identifier
213 * \param[in] box The box tensor
215 void set_pbc(t_pbc* pbc, PbcType pbcType, const matrix box);
217 /*! \brief Initiate the periodic boundary condition algorithms.
219 * As set_pbc, but additionally sets that correct distances can
220 * be obtained using (combinations of) single box-vector shifts.
221 * Should be used with pbc_dx_aiuc.
222 * If domdecCells!=NULL pbc is not used for directions
223 * with dd->nc[i]==1 with bSingleDir==TRUE or
224 * with dd->nc[i]<=2 with bSingleDir==FALSE.
225 * Note that when no PBC is required only pbc->pbcType is set,
226 * the rest of the struct will be invalid.
228 * \param[in,out] pbc The pbc information structure
229 * \param[in] pbcType The PBC identifier
230 * \param[in] domdecCells 3D integer vector describing the number of DD cells
231 * or nullptr if not using DD.
232 * \param[in] bSingleDir TRUE if DD communicates only in one direction along dimensions
233 * \param[in] box The box tensor
234 * \return the pbc structure when pbc operations are required, NULL otherwise.
236 t_pbc* set_pbc_dd(t_pbc* pbc, PbcType pbcType, const ivec domdecCells, bool bSingleDir, const matrix box);
238 /*! \brief Compute distance with PBC
240 * Calculate the correct distance vector from x2 to x1 and put it in dx.
241 * set_pbc must be called before ever calling this routine.
243 * Note that for triclinic boxes that do not obey the GROMACS unit-cell
244 * restrictions, pbc_dx and pbc_dx_aiuc will not correct for PBC.
245 * \param[in,out] pbc The pbc information structure
246 * \param[in] x1 Coordinates for particle 1
247 * \param[in] x2 Coordinates for particle 2
248 * \param[out] dx Distance vector
250 void pbc_dx(const t_pbc* pbc, const rvec x1, const rvec x2, rvec dx);
252 /*! \brief Compute distance vector for simple PBC types
254 * Calculate the correct distance vector from x2 to x1 and put it in dx,
255 * This function can only be used when all atoms are in the rectangular
256 * or triclinic unit-cell.
257 * set_pbc_dd or set_pbc must be called before ever calling this routine.
258 * \param[in,out] pbc The pbc information structure
259 * \param[in] x1 Coordinates for particle 1
260 * \param[in] x2 Coordinates for particle 2
261 * \param[out] dx Distance vector
262 * \return the ishift required to shift x1 at closest distance to x2;
263 * i.e. if 0<=ishift<SHIFTS then x1 - x2 + shift_vec[ishift] = dx
264 * (see calc_shifts below on how to obtain shift_vec)
266 int pbc_dx_aiuc(const t_pbc* pbc, const rvec x1, const rvec x2, rvec dx);
268 /*! \brief Compute distance with PBC
270 * As pbc_dx, but for double precision vectors.
271 * set_pbc must be called before ever calling this routine.
272 * \param[in,out] pbc The pbc information structure
273 * \param[in] x1 Coordinates for particle 1
274 * \param[in] x2 Coordinates for particle 2
275 * \param[out] dx Distance vector
277 void pbc_dx_d(const t_pbc* pbc, const dvec x1, const dvec x2, dvec dx);
279 /*! \brief Computes shift vectors
281 * This routine calculates ths shift vectors necessary to use the
282 * neighbor searching routine.
283 * \param[in] box The simulation box
284 * \param[out] shift_vec The shifting vectors
286 void calc_shifts(const matrix box, gmx::ArrayRef<gmx::RVec> shift_vec);
288 /*! \brief Calculates the center of the box.
290 * See the description for the enum ecenter above.
291 * \param[in] ecenter Description of center type
292 * \param[in] box The simulation box
293 * \param[out] box_center The center of the box
295 void calc_box_center(int ecenter, const matrix box, rvec box_center);
297 /*! \brief Calculates the NTRICIMG box images
299 * \param[in] box The simulation box
300 * \param[out] img The triclinic box images
302 void calc_triclinic_images(const matrix box, rvec img[]);
304 /*! \brief Calculates the NCUCVERT vertices of a compact unitcell
306 * \param[in] ecenter The center type
307 * \param[in] box The simulation box
308 * \param[out] vert The vertices
310 void calc_compact_unitcell_vertices(int ecenter, const matrix box, rvec vert[]);
312 /*! \brief Compute unitcell edges
314 * \return an array of unitcell edges of length NCUCEDGE*2,
315 * this is an index in vert[], which is calculated by calc_unitcell_vertices.
316 * The index consists of NCUCEDGE pairs of vertex indices.
317 * The index does not change, so it needs to be retrieved only once.
319 int* compact_unitcell_edges();
321 /*! \brief Put atoms inside the simulations box
323 * These routines puts ONE or ALL atoms in the box, not caring
324 * about charge groups!
325 * Also works for triclinic cells.
327 * \param[in] pbcType The pbc type
328 * \param[in] box The simulation box
329 * \param[in,out] x The coordinates of the atoms
331 void put_atoms_in_box(PbcType pbcType, const matrix box, gmx::ArrayRef<gmx::RVec> x);
333 /*! \brief Parallellizes put_atoms_in_box()
335 * This wrapper function around put_atoms_in_box() with the ugly manual
336 * workload splitting is needed to avoid silently introducing multithreading
339 * \param[in] pbcType The pbc type
340 * \param[in] box The simulation box
341 * \param[in,out] x The coordinates of the atoms
342 * \param[in] nth number of threads to be used in the given module
344 void put_atoms_in_box_omp(PbcType pbcType, const matrix box, gmx::ArrayRef<gmx::RVec> x, gmx_unused int nth);
346 /*! \brief Put atoms inside triclinic box
348 * This puts ALL atoms in the triclinic unit cell, centered around the
349 * box center as calculated by calc_box_center.
350 * \param[in] ecenter The pbc center type
351 * \param[in] box The simulation box
352 * \param[in,out] x The coordinates of the atoms
354 void put_atoms_in_triclinic_unitcell(int ecenter, const matrix box, gmx::ArrayRef<gmx::RVec> x);
356 /*! \brief Put atoms inside the unitcell
358 * This puts ALL atoms at the closest distance for the center of the box
359 * as calculated by calc_box_center.
360 * When \p pbcType=PbcType::Unset, the type of pbc is guessed from the box matrix.
362 * \param[in] pbcType The pbc type
363 * \param[in] ecenter The pbc center type
364 * \param[in] box The simulation box
365 * \param[in,out] x The coordinates of the atoms
367 void put_atoms_in_compact_unitcell(PbcType pbcType, int ecenter, const matrix box, gmx::ArrayRef<gmx::RVec> x);
369 /*! \brief Make all molecules whole by shifting positions
371 * \param[in] fplog Log file
372 * \param[in] pbcType The PBC type
373 * \param[in] box The simulation box
374 * \param[in] mtop System topology definition
375 * \param[in,out] x The coordinates of the atoms
377 void do_pbc_first_mtop(FILE* fplog, PbcType pbcType, const matrix box, const gmx_mtop_t* mtop, rvec x[]);
379 /*! \brief Make molecules consisting of multiple charge groups whole by shifting positions
381 * \param[in] pbcType The PBC type
382 * \param[in] box The simulation box
383 * \param[in] mtop System topology definition
384 * \param[in,out] x The coordinates of the atoms
386 void do_pbc_mtop(PbcType pbcType, const matrix box, const gmx_mtop_t* mtop, rvec x[]);