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
135 /*! \brief Returns the number of dimensions that use pbc
137 * \param[in] pbcType The periodic boundary condition type
138 * \return the number of dimensions that use pbc, starting at X
140 int numPbcDimensions(PbcType pbcType);
142 /*! \brief Dump the contents of the pbc structure to the file
144 * \param[in] fp The file pointer to write to
145 * \param[in] pbc The periodic boundary condition information structure
147 void dump_pbc(FILE* fp, t_pbc* pbc);
149 /*! \brief Check the box for consistency
151 * When \p pbcType=PbcTypes::Unset, the type of pbc is guessed from the box matrix.
153 * \param[in] pbcType The pbc identifier
154 * \param[in] box The box matrix
155 * \return NULL if the box is supported by Gromacs.
156 * Otherwise returns a string with the problem.
158 const char* check_box(PbcType pbcType, const matrix box);
160 /*! \brief Creates box matrix from edge lengths and angles.
162 * \param[in,out] box The box matrix
163 * \param[in] vec The edge lengths
164 * \param[in] angleInDegrees The angles
166 void matrix_convert(matrix box, const rvec vec, const rvec angleInDegrees);
168 /*! \brief Compute the maximum cutoff for the box
170 * Returns the square of the maximum cut-off allowed for the box,
171 * taking into account that the grid neighborsearch code and pbc_dx
172 * only check combinations of single box-vector shifts.
174 * \param[in] pbcType The pbc identifier
175 * \param[in] box The box matrix
176 * \return the maximum cut-off.
178 real max_cutoff2(PbcType pbcType, const matrix box);
180 /*! \brief Guess PBC type
182 * Guesses the type of periodic boundary conditions using the box
184 * \param[in] box The box matrix
185 * \return The pbc type identifier
187 PbcType guessPbcType(const matrix box);
189 /*! \brief Corrects the box if necessary
191 * Checks for un-allowed box angles and corrects the box.
193 * \param[in] fplog File for debug output
194 * \param[in] step The MD step number
195 * \param[in] box The simulation cell
196 * \return TRUE when the box was corrected.
198 gmx_bool correct_box(FILE* fplog, int step, tensor box);
200 /*! \brief Initiate the periodic boundary condition algorithms.
202 * pbc_dx will not use pbc and return the normal difference vector
203 * when one or more of the diagonal elements of box are zero.
204 * When \p pbcType=PbcType::Unset, the type of pbc is guessed from the box matrix.
206 * \param[in,out] pbc The pbc information structure
207 * \param[in] pbcType The PBC identifier
208 * \param[in] box The box tensor
210 void set_pbc(t_pbc* pbc, PbcType pbcType, const matrix box);
212 /*! \brief Initiate the periodic boundary condition algorithms.
214 * As set_pbc, but additionally sets that correct distances can
215 * be obtained using (combinations of) single box-vector shifts.
216 * Should be used with pbc_dx_aiuc.
217 * If domdecCells!=NULL pbc is not used for directions
218 * with dd->nc[i]==1 with bSingleDir==TRUE or
219 * with dd->nc[i]<=2 with bSingleDir==FALSE.
220 * Note that when no PBC is required only pbc->pbcType is set,
221 * the rest of the struct will be invalid.
223 * \param[in,out] pbc The pbc information structure
224 * \param[in] pbcType The PBC identifier
225 * \param[in] domdecCells 3D integer vector describing the number of DD cells
226 * or nullptr if not using DD.
227 * \param[in] bSingleDir TRUE if DD communicates only in one direction along dimensions
228 * \param[in] box The box tensor
229 * \return the pbc structure when pbc operations are required, NULL otherwise.
231 t_pbc* set_pbc_dd(t_pbc* pbc, PbcType pbcType, const ivec domdecCells, gmx_bool bSingleDir, const matrix box);
233 /*! \brief Compute distance with PBC
235 * Calculate the correct distance vector from x2 to x1 and put it in dx.
236 * set_pbc must be called before ever calling this routine.
238 * Note that for triclinic boxes that do not obey the GROMACS unit-cell
239 * restrictions, pbc_dx and pbc_dx_aiuc will not correct for PBC.
240 * \param[in,out] pbc The pbc information structure
241 * \param[in] x1 Coordinates for particle 1
242 * \param[in] x2 Coordinates for particle 2
243 * \param[out] dx Distance vector
245 void pbc_dx(const t_pbc* pbc, const rvec x1, const rvec x2, rvec dx);
247 /*! \brief Compute distance vector for simple PBC types
249 * Calculate the correct distance vector from x2 to x1 and put it in dx,
250 * This function can only be used when all atoms are in the rectangular
251 * or triclinic unit-cell.
252 * set_pbc_dd or set_pbc must be called before ever calling this routine.
253 * \param[in,out] pbc The pbc information structure
254 * \param[in] x1 Coordinates for particle 1
255 * \param[in] x2 Coordinates for particle 2
256 * \param[out] dx Distance vector
257 * \return the ishift required to shift x1 at closest distance to x2;
258 * i.e. if 0<=ishift<SHIFTS then x1 - x2 + shift_vec[ishift] = dx
259 * (see calc_shifts below on how to obtain shift_vec)
261 int pbc_dx_aiuc(const t_pbc* pbc, const rvec x1, const rvec x2, rvec dx);
263 /*! \brief Compute distance with PBC
265 * As pbc_dx, but for double precision vectors.
266 * set_pbc must be called before ever calling this routine.
267 * \param[in,out] pbc The pbc information structure
268 * \param[in] x1 Coordinates for particle 1
269 * \param[in] x2 Coordinates for particle 2
270 * \param[out] dx Distance vector
272 void pbc_dx_d(const t_pbc* pbc, const dvec x1, const dvec x2, dvec dx);
274 /*! \brief Computes shift vectors
276 * This routine calculates ths shift vectors necessary to use the
277 * neighbor searching routine.
278 * \param[in] box The simulation box
279 * \param[out] shift_vec The shifting vectors
281 void calc_shifts(const matrix box, rvec shift_vec[]);
283 /*! \brief Calculates the center of the box.
285 * See the description for the enum ecenter above.
286 * \param[in] ecenter Description of center type
287 * \param[in] box The simulation box
288 * \param[out] box_center The center of the box
290 void calc_box_center(int ecenter, const matrix box, rvec box_center);
292 /*! \brief Calculates the NTRICIMG box images
294 * \param[in] box The simulation box
295 * \param[out] img The triclinic box images
297 void calc_triclinic_images(const matrix box, rvec img[]);
299 /*! \brief Calculates the NCUCVERT vertices of a compact unitcell
301 * \param[in] ecenter The center type
302 * \param[in] box The simulation box
303 * \param[out] vert The vertices
305 void calc_compact_unitcell_vertices(int ecenter, const matrix box, rvec vert[]);
307 /*! \brief Compute unitcell edges
309 * \return an array of unitcell edges of length NCUCEDGE*2,
310 * this is an index in vert[], which is calculated by calc_unitcell_vertices.
311 * The index consists of NCUCEDGE pairs of vertex indices.
312 * The index does not change, so it needs to be retrieved only once.
314 int* compact_unitcell_edges();
316 /*! \brief Put atoms inside the simulations box
318 * These routines puts ONE or ALL atoms in the box, not caring
319 * about charge groups!
320 * Also works for triclinic cells.
322 * \param[in] pbcType The pbc type
323 * \param[in] box The simulation box
324 * \param[in,out] x The coordinates of the atoms
326 void put_atoms_in_box(PbcType pbcType, const matrix box, gmx::ArrayRef<gmx::RVec> x);
328 /*! \brief Parallellizes put_atoms_in_box()
330 * This wrapper function around put_atoms_in_box() with the ugly manual
331 * workload splitting is needed to avoid silently introducing multithreading
334 * \param[in] pbcType The pbc type
335 * \param[in] box The simulation box
336 * \param[in,out] x The coordinates of the atoms
337 * \param[in] nth number of threads to be used in the given module
339 void put_atoms_in_box_omp(PbcType pbcType, const matrix box, gmx::ArrayRef<gmx::RVec> x, gmx_unused int nth);
341 /*! \brief Put atoms inside triclinic box
343 * This puts ALL atoms in the triclinic unit cell, centered around the
344 * box center as calculated by calc_box_center.
345 * \param[in] ecenter The pbc center type
346 * \param[in] box The simulation box
347 * \param[in,out] x The coordinates of the atoms
349 void put_atoms_in_triclinic_unitcell(int ecenter, const matrix box, gmx::ArrayRef<gmx::RVec> x);
351 /*! \brief Put atoms inside the unitcell
353 * This puts ALL atoms at the closest distance for the center of the box
354 * as calculated by calc_box_center.
355 * When \p pbcType=PbcType::Unset, the type of pbc is guessed from the box matrix.
357 * \param[in] pbcType The pbc type
358 * \param[in] ecenter The pbc center type
359 * \param[in] box The simulation box
360 * \param[in,out] x The coordinates of the atoms
362 void put_atoms_in_compact_unitcell(PbcType pbcType, int ecenter, const matrix box, gmx::ArrayRef<gmx::RVec> x);
364 /*! \brief Make all molecules whole by shifting positions
366 * \param[in] fplog Log file
367 * \param[in] pbcType The PBC type
368 * \param[in] box The simulation box
369 * \param[in] mtop System topology definition
370 * \param[in,out] x The coordinates of the atoms
372 void do_pbc_first_mtop(FILE* fplog, PbcType pbcType, const matrix box, const gmx_mtop_t* mtop, rvec x[]);
374 /*! \brief Make molecules consisting of multiple charge groups whole by shifting positions
376 * \param[in] pbcType The PBC type
377 * \param[in] box The simulation box
378 * \param[in] mtop System topology definition
379 * \param[in,out] x The coordinates of the atoms
381 void do_pbc_mtop(PbcType pbcType, const matrix box, const gmx_mtop_t* mtop, rvec x[]);