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, 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/arrayref.h"
47 #include "gromacs/utility/basedefinitions.h"
48 #include "gromacs/utility/enumerationhelpers.h"
49 #include "gromacs/utility/real.h"
54 //! Enumeration that contains all supported periodic boundary setups.
55 enum class PbcType : int
57 Xyz = 0, //!< Periodic boundaries in all dimensions.
58 No = 1, //!< No periodic boundaries.
59 XY = 2, //!< Only two dimensions are periodic.
60 Screw = 3, //!< Screw.
61 Unset = 4, //!< The type of PBC is not set or invalid.
65 //! Names for all values in PBC types enumeration
66 extern const gmx::EnumerationArray<PbcType, std::string> c_pbcTypeNames;
68 /* Maximum number of combinations of single triclinic box vectors
69 * required to shift atoms that are within a brick of the size of
70 * the diagonal of the box to within the maximum cut-off distance.
72 #define MAX_NTRICVEC 12
74 /*! \brief Structure containing info on periodic boundary conditions */
79 //! Number of dimensions in which PBC is exerted
81 /*! \brief Determines how to compute distance vectors.
83 * Indicator of how to compute distance vectors, depending
84 * on PBC type (depends on pbcType and dimensions with(out) DD)
88 /*! \brief Used for selecting which dimensions to use in PBC.
90 * In case of 1-D PBC this indicates which dimension is used,
91 * in case of 2-D PBC this indicates the opposite
94 //! The simulation box
96 //! The lengths of the diagonal of the full box
98 //! Halve of the above
100 //! Negative of the above
102 //! Maximum allowed cutoff squared for the box and PBC used
104 /*! \brief Number of triclinic shift vectors.
106 * Number of triclinic shift vectors depends on the skewedness
107 * of the box, that is mostly on the angles. For triclinic boxes
108 * we first take the closest image along each Cartesian dimension
109 * independently. When the resulting distance^2 is larger than
110 * max_cutoff2, up to ntric_vec triclinic shift vectors need to
111 * be tried. Because of the restrictions imposed on the unit-cell
112 * by GROMACS, ntric_vec <= MAX_NTRICVEC = 12.
115 //! The triclinic shift vectors in grid cells. Internal use only.
116 ivec tric_shift[MAX_NTRICVEC];
117 //! The triclinic shift vectors in length units
118 rvec tric_vec[MAX_NTRICVEC];
121 #define TRICLINIC(box) ((box)[YY][XX] != 0 || (box)[ZZ][XX] != 0 || (box)[ZZ][YY] != 0)
129 ecenterTRIC, /* 0.5*(a+b+c) */
130 ecenterRECT, /* (0.5*a[x],0.5*b[y],0.5*c[z]) */
131 ecenterZERO, /* (0,0,0) */
132 ecenterDEF = ecenterTRIC
137 /*! \brief Returns the number of dimensions that use pbc
139 * \param[in] pbcType The periodic boundary condition type
140 * \return the number of dimensions that use pbc, starting at X
142 int numPbcDimensions(PbcType pbcType);
144 /*! \brief Dump the contents of the pbc structure to the file
146 * \param[in] fp The file pointer to write to
147 * \param[in] pbc The periodic boundary condition information structure
149 void dump_pbc(FILE* fp, t_pbc* pbc);
151 /*! \brief Check the box for consistency
153 * When \p pbcType=PbcTypes::Unset, the type of pbc is guessed from the box matrix.
155 * \param[in] pbcType The pbc identifier
156 * \param[in] box The box matrix
157 * \return NULL if the box is supported by Gromacs.
158 * Otherwise returns a string with the problem.
160 const char* check_box(PbcType pbcType, const matrix box);
162 /*! \brief Creates box matrix from edge lengths and angles.
164 * \param[in,out] box The box matrix
165 * \param[in] vec The edge lengths
166 * \param[in] angleInDegrees The angles
168 void matrix_convert(matrix box, const rvec vec, const rvec angleInDegrees);
170 /*! \brief Compute the maximum cutoff for the box
172 * Returns the square of the maximum cut-off allowed for the box,
173 * taking into account that the grid neighborsearch code and pbc_dx
174 * only check combinations of single box-vector shifts.
176 * \param[in] pbcType The pbc identifier
177 * \param[in] box The box matrix
178 * \return the maximum cut-off.
180 real max_cutoff2(PbcType pbcType, const matrix box);
182 /*! \brief Guess PBC type
184 * Guesses the type of periodic boundary conditions using the box
186 * \param[in] box The box matrix
187 * \return The pbc type identifier
189 PbcType guessPbcType(const matrix box);
191 /*! \brief Corrects the box if necessary
193 * Checks for un-allowed box angles and corrects the box
194 * and the integer shift vectors in the graph (if \p graph!=NULL) if necessary.
196 * \param[in] fplog File for debug output
197 * \param[in] step The MD step number
198 * \param[in] box The simulation cell
199 * \param[in] graph Information about molecular connectivity
200 * \return TRUE when the box was corrected.
202 gmx_bool correct_box(FILE* fplog, int step, tensor box, struct t_graph* graph);
204 /*! \brief Initiate the periodic boundary condition algorithms.
206 * pbc_dx will not use pbc and return the normal difference vector
207 * when one or more of the diagonal elements of box are zero.
208 * When \p pbcType=PbcType::Unset, the type of pbc is guessed from the box matrix.
210 * \param[in,out] pbc The pbc information structure
211 * \param[in] pbcType The PBC identifier
212 * \param[in] box The box tensor
214 void set_pbc(t_pbc* pbc, PbcType pbcType, const matrix box);
216 /*! \brief Initiate the periodic boundary condition algorithms.
218 * As set_pbc, but additionally sets that correct distances can
219 * be obtained using (combinations of) single box-vector shifts.
220 * Should be used with pbc_dx_aiuc.
221 * If domdecCells!=NULL pbc is not used for directions
222 * with dd->nc[i]==1 with bSingleDir==TRUE or
223 * with dd->nc[i]<=2 with bSingleDir==FALSE.
224 * Note that when no PBC is required only pbc->pbcType is set,
225 * the rest of the struct will be invalid.
227 * \param[in,out] pbc The pbc information structure
228 * \param[in] pbcType The PBC identifier
229 * \param[in] domdecCells 3D integer vector describing the number of DD cells
230 * or nullptr if not using DD.
231 * \param[in] bSingleDir TRUE if DD communicates only in one direction along dimensions
232 * \param[in] box The box tensor
233 * \return the pbc structure when pbc operations are required, NULL otherwise.
235 t_pbc* set_pbc_dd(t_pbc* pbc, PbcType pbcType, const ivec domdecCells, gmx_bool bSingleDir, const matrix box);
237 /*! \brief Compute distance with PBC
239 * Calculate the correct distance vector from x2 to x1 and put it in dx.
240 * set_pbc must be called before ever calling this routine.
242 * Note that for triclinic boxes that do not obey the GROMACS unit-cell
243 * restrictions, pbc_dx and pbc_dx_aiuc will not correct for PBC.
244 * \param[in,out] pbc The pbc information structure
245 * \param[in] x1 Coordinates for particle 1
246 * \param[in] x2 Coordinates for particle 2
247 * \param[out] dx Distance vector
249 void pbc_dx(const t_pbc* pbc, const rvec x1, const rvec x2, rvec dx);
251 /*! \brief Compute distance vector for simple PBC types
253 * Calculate the correct distance vector from x2 to x1 and put it in dx,
254 * This function can only be used when all atoms are in the rectangular
255 * or triclinic unit-cell.
256 * set_pbc_dd or set_pbc must be called before ever calling this routine.
257 * \param[in,out] pbc The pbc information structure
258 * \param[in] x1 Coordinates for particle 1
259 * \param[in] x2 Coordinates for particle 2
260 * \param[out] dx Distance vector
261 * \return the ishift required to shift x1 at closest distance to x2;
262 * i.e. if 0<=ishift<SHIFTS then x1 - x2 + shift_vec[ishift] = dx
263 * (see calc_shifts below on how to obtain shift_vec)
265 int pbc_dx_aiuc(const t_pbc* pbc, const rvec x1, const rvec x2, rvec dx);
267 /*! \brief Compute distance with PBC
269 * As pbc_dx, but for double precision vectors.
270 * set_pbc must be called before ever calling this routine.
271 * \param[in,out] pbc The pbc information structure
272 * \param[in] x1 Coordinates for particle 1
273 * \param[in] x2 Coordinates for particle 2
274 * \param[out] dx Distance vector
276 void pbc_dx_d(const t_pbc* pbc, const dvec x1, const dvec x2, dvec dx);
278 /*! \brief Computes shift vectors
280 * This routine calculates ths shift vectors necessary to use the
281 * neighbor searching routine.
282 * \param[in] box The simulation box
283 * \param[out] shift_vec The shifting vectors
285 void calc_shifts(const matrix box, rvec shift_vec[]);
287 /*! \brief Calculates the center of the box.
289 * See the description for the enum ecenter above.
290 * \param[in] ecenter Description of center type
291 * \param[in] box The simulation box
292 * \param[out] box_center The center of the box
294 void calc_box_center(int ecenter, const matrix box, rvec box_center);
296 /*! \brief Calculates the NTRICIMG box images
298 * \param[in] box The simulation box
299 * \param[out] img The triclinic box images
301 void calc_triclinic_images(const matrix box, rvec img[]);
303 /*! \brief Calculates the NCUCVERT vertices of a compact unitcell
305 * \param[in] ecenter The center type
306 * \param[in] box The simulation box
307 * \param[out] vert The vertices
309 void calc_compact_unitcell_vertices(int ecenter, const matrix box, rvec vert[]);
311 /*! \brief Compute unitcell edges
313 * \return an array of unitcell edges of length NCUCEDGE*2,
314 * this is an index in vert[], which is calculated by calc_unitcell_vertices.
315 * The index consists of NCUCEDGE pairs of vertex indices.
316 * The index does not change, so it needs to be retrieved only once.
318 int* compact_unitcell_edges();
320 /*! \brief Put atoms inside the simulations box
322 * These routines puts ONE or ALL atoms in the box, not caring
323 * about charge groups!
324 * Also works for triclinic cells.
326 * \param[in] pbcType The pbc type
327 * \param[in] box The simulation box
328 * \param[in,out] x The coordinates of the atoms
330 void put_atoms_in_box(PbcType pbcType, const matrix box, gmx::ArrayRef<gmx::RVec> x);
332 /*! \brief Parallellizes put_atoms_in_box()
334 * This wrapper function around put_atoms_in_box() with the ugly manual
335 * workload splitting is needed to avoid silently introducing multithreading
338 * \param[in] pbcType The pbc type
339 * \param[in] box The simulation box
340 * \param[in,out] x The coordinates of the atoms
341 * \param[in] nth number of threads to be used in the given module
343 void put_atoms_in_box_omp(PbcType pbcType, const matrix box, gmx::ArrayRef<gmx::RVec> x, gmx_unused int nth);
345 /*! \brief Put atoms inside triclinic box
347 * This puts ALL atoms in the triclinic unit cell, centered around the
348 * box center as calculated by calc_box_center.
349 * \param[in] ecenter The pbc center type
350 * \param[in] box The simulation box
351 * \param[in,out] x The coordinates of the atoms
353 void put_atoms_in_triclinic_unitcell(int ecenter, const matrix box, gmx::ArrayRef<gmx::RVec> x);
355 /*! \brief Put atoms inside the unitcell
357 * This puts ALL atoms at the closest distance for the center of the box
358 * as calculated by calc_box_center.
359 * When \p pbcType=PbcType::Unset, the type of pbc is guessed from the box matrix.
361 * \param[in] pbcType The pbc type
362 * \param[in] ecenter The pbc center type
363 * \param[in] box The simulation box
364 * \param[in,out] x The coordinates of the atoms
366 void put_atoms_in_compact_unitcell(PbcType pbcType, int ecenter, const matrix box, gmx::ArrayRef<gmx::RVec> x);
368 /*! \brief Make all molecules whole by shifting positions
370 * \param[in] fplog Log file
371 * \param[in] pbcType The PBC type
372 * \param[in] box The simulation box
373 * \param[in] mtop System topology definition
374 * \param[in,out] x The coordinates of the atoms
376 void do_pbc_first_mtop(FILE* fplog, PbcType pbcType, const matrix box, const gmx_mtop_t* mtop, rvec x[]);
378 /*! \brief Make molecules consisting of multiple charge groups whole by shifting positions
380 * \param[in] pbcType The PBC type
381 * \param[in] box The simulation box
382 * \param[in] mtop System topology definition
383 * \param[in,out] x The coordinates of the atoms
385 void do_pbc_mtop(PbcType pbcType, const matrix box, const gmx_mtop_t* mtop, rvec x[]);