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
43 #include "gromacs/math/vectypes.h"
44 #include "gromacs/utility/arrayref.h"
45 #include "gromacs/utility/basedefinitions.h"
46 #include "gromacs/utility/real.h"
60 //! Strings corresponding to epbc enum values.
61 extern const char* epbc_names[epbcNR + 1];
63 /* Maximum number of combinations of single triclinic box vectors
64 * required to shift atoms that are within a brick of the size of
65 * the diagonal of the box to within the maximum cut-off distance.
67 #define MAX_NTRICVEC 12
69 /*! \brief Structure containing info on periodic boundary conditions */
74 //! Number of dimensions in which PBC is exerted
76 /*! \brief Determines how to compute distance vectors.
78 * Indicator of how to compute distance vectors, depending
79 * on PBC type (depends on ePBC and dimensions with(out) DD)
83 /*! \brief Used for selecting which dimensions to use in PBC.
85 * In case of 1-D PBC this indicates which dimension is used,
86 * in case of 2-D PBC this indicates the opposite
89 //! The simulation box
91 //! The lengths of the diagonal of the full box
93 //! Halve of the above
95 //! Negative of the above
97 //! Maximum allowed cutoff squared for the box and PBC used
99 /*! \brief Number of triclinic shift vectors.
101 * Number of triclinic shift vectors depends on the skewedness
102 * of the box, that is mostly on the angles. For triclinic boxes
103 * we first take the closest image along each Cartesian dimension
104 * independently. When the resulting distance^2 is larger than
105 * max_cutoff2, up to ntric_vec triclinic shift vectors need to
106 * be tried. Because of the restrictions imposed on the unit-cell
107 * by GROMACS, ntric_vec <= MAX_NTRICVEC = 12.
110 //! The triclinic shift vectors in grid cells. Internal use only.
111 ivec tric_shift[MAX_NTRICVEC];
112 //! The triclinic shift vectors in length units
113 rvec tric_vec[MAX_NTRICVEC];
116 #define TRICLINIC(box) ((box)[YY][XX] != 0 || (box)[ZZ][XX] != 0 || (box)[ZZ][YY] != 0)
124 ecenterTRIC, /* 0.5*(a+b+c) */
125 ecenterRECT, /* (0.5*a[x],0.5*b[y],0.5*c[z]) */
126 ecenterZERO, /* (0,0,0) */
127 ecenterDEF = ecenterTRIC
132 /*! \brief Returns the number of dimensions that use pbc
134 * \param[in] ePBC The periodic boundary condition type
135 * \return the number of dimensions that use pbc, starting at X
137 int ePBC2npbcdim(int ePBC);
139 /*! \brief Dump the contents of the pbc structure to the file
141 * \param[in] fp The file pointer to write to
142 * \param[in] pbc The periodic boundary condition information structure
144 void dump_pbc(FILE* fp, t_pbc* pbc);
146 /*! \brief Check the box for consistency
148 * \param[in] ePBC The pbc identifier
149 * \param[in] box The box matrix
150 * \return NULL if the box is supported by Gromacs.
151 * Otherwise returns a string with the problem.
152 * When ePBC=-1, the type of pbc is guessed from the box matrix.
154 const char* check_box(int ePBC, const matrix box);
156 /*! \brief Creates box matrix from edge lengths and angles.
158 * \param[in,out] box The box matrix
159 * \param[in] vec The edge lengths
160 * \param[in] angleInDegrees The angles
162 void matrix_convert(matrix box, const rvec vec, const rvec angleInDegrees);
164 /*! \brief Compute the maximum cutoff for the box
166 * Returns the square of the maximum cut-off allowed for the box,
167 * taking into account that the grid neighborsearch code and pbc_dx
168 * only check combinations of single box-vector shifts.
169 * \param[in] ePBC The pbc identifier
170 * \param[in] box The box matrix
171 * \return the maximum cut-off.
173 real max_cutoff2(int ePBC, const matrix box);
175 /*! \brief Guess PBC typr
177 * Guesses the type of periodic boundary conditions using the box
178 * \param[in] box The box matrix
179 * \return The pbc identifier
181 int guess_ePBC(const matrix box);
183 /*! \brief Corrects the box if necessary
185 * Checks for un-allowed box angles and corrects the box
186 * and the integer shift vectors in the graph (if graph!=NULL) if necessary.
187 * \param[in] fplog File for debug output
188 * \param[in] step The MD step number
189 * \param[in] box The simulation cell
190 * \param[in] graph Information about molecular connectivity
191 * \return TRUE when the box was corrected.
193 gmx_bool correct_box(FILE* fplog, int step, tensor box, struct t_graph* graph);
195 /*! \brief Initiate the periodic boundary condition algorithms.
197 * pbc_dx will not use pbc and return the normal difference vector
198 * when one or more of the diagonal elements of box are zero.
199 * When ePBC=-1, the type of pbc is guessed from the box matrix.
200 * \param[in,out] pbc The pbc information structure
201 * \param[in] ePBC The PBC identifier
202 * \param[in] box The box tensor
204 void set_pbc(t_pbc* pbc, int ePBC, const matrix box);
206 /*! \brief Initiate the periodic boundary condition algorithms.
208 * As set_pbc, but additionally sets that correct distances can
209 * be obtained using (combinations of) single box-vector shifts.
210 * Should be used with pbc_dx_aiuc.
211 * If domdecCells!=NULL pbc is not used for directions
212 * with dd->nc[i]==1 with bSingleDir==TRUE or
213 * with dd->nc[i]<=2 with bSingleDir==FALSE.
214 * Note that when no PBC is required only pbc->ePBC is set,
215 * the rest of the struct will be invalid.
216 * \param[in,out] pbc The pbc information structure
217 * \param[in] ePBC The PBC identifier
218 * \param[in] domdecCells 3D integer vector describing the number of DD cells
219 * or nullptr if not using DD.
220 * \param[in] bSingleDir TRUE if DD communicates only in one direction along dimensions
221 * \param[in] box The box tensor
222 * \return the pbc structure when pbc operations are required, NULL otherwise.
224 t_pbc* set_pbc_dd(t_pbc* pbc, int ePBC, const ivec domdecCells, gmx_bool bSingleDir, const matrix box);
226 /*! \brief Compute distance with PBC
228 * Calculate the correct distance vector from x2 to x1 and put it in dx.
229 * set_pbc must be called before ever calling this routine.
231 * Note that for triclinic boxes that do not obey the GROMACS unit-cell
232 * restrictions, pbc_dx and pbc_dx_aiuc will not correct for PBC.
233 * \param[in,out] pbc The pbc information structure
234 * \param[in] x1 Coordinates for particle 1
235 * \param[in] x2 Coordinates for particle 2
236 * \param[out] dx Distance vector
238 void pbc_dx(const t_pbc* pbc, const rvec x1, const rvec x2, rvec dx);
240 /*! \brief Compute distance vector for simple PBC types
242 * Calculate the correct distance vector from x2 to x1 and put it in dx,
243 * This function can only be used when all atoms are in the rectangular
244 * or triclinic unit-cell.
245 * set_pbc_dd or set_pbc must be called before ever calling this routine.
246 * \param[in,out] pbc The pbc information structure
247 * \param[in] x1 Coordinates for particle 1
248 * \param[in] x2 Coordinates for particle 2
249 * \param[out] dx Distance vector
250 * \return the ishift required to shift x1 at closest distance to x2;
251 * i.e. if 0<=ishift<SHIFTS then x1 - x2 + shift_vec[ishift] = dx
252 * (see calc_shifts below on how to obtain shift_vec)
254 int pbc_dx_aiuc(const t_pbc* pbc, const rvec x1, const rvec x2, rvec dx);
256 /*! \brief Compute distance with PBC
258 * As pbc_dx, but for double precision vectors.
259 * set_pbc must be called before ever calling this routine.
260 * \param[in,out] pbc The pbc information structure
261 * \param[in] x1 Coordinates for particle 1
262 * \param[in] x2 Coordinates for particle 2
263 * \param[out] dx Distance vector
265 void pbc_dx_d(const t_pbc* pbc, const dvec x1, const dvec x2, dvec dx);
267 /*! \brief Computes shift vectors
269 * This routine calculates ths shift vectors necessary to use the
270 * neighbor searching routine.
271 * \param[in] box The simulation box
272 * \param[out] shift_vec The shifting vectors
274 void calc_shifts(const matrix box, rvec shift_vec[]);
276 /*! \brief Calculates the center of the box.
278 * See the description for the enum ecenter above.
279 * \param[in] ecenter Description of center type
280 * \param[in] box The simulation box
281 * \param[out] box_center The center of the box
283 void calc_box_center(int ecenter, const matrix box, rvec box_center);
285 /*! \brief Calculates the NTRICIMG box images
287 * \param[in] box The simulation box
288 * \param[out] img The triclinic box images
290 void calc_triclinic_images(const matrix box, rvec img[]);
292 /*! \brief Calculates the NCUCVERT vertices of a compact unitcell
294 * \param[in] ecenter The center type
295 * \param[in] box The simulation box
296 * \param[out] vert The vertices
298 void calc_compact_unitcell_vertices(int ecenter, const matrix box, rvec vert[]);
300 /*! \brief Compute unitcell edges
302 * \return an array of unitcell edges of length NCUCEDGE*2,
303 * this is an index in vert[], which is calculated by calc_unitcell_vertices.
304 * The index consists of NCUCEDGE pairs of vertex indices.
305 * The index does not change, so it needs to be retrieved only once.
307 int* compact_unitcell_edges();
309 /*! \brief Put atoms inside the simulations box
311 * These routines puts ONE or ALL atoms in the box, not caring
312 * about charge groups!
313 * Also works for triclinic cells.
314 * \param[in] ePBC The pbc type
315 * \param[in] box The simulation box
316 * \param[in,out] x The coordinates of the atoms
318 void put_atoms_in_box(int ePBC, const matrix box, gmx::ArrayRef<gmx::RVec> x);
320 /*! \brief Parallellizes put_atoms_in_box()
322 * This wrapper function around put_atoms_in_box() with the ugly manual
323 * workload splitting is needed to avoid silently introducing multithreading
325 * \param[in] ePBC The pbc type
326 * \param[in] box The simulation box
327 * \param[in,out] x The coordinates of the atoms
328 * \param[in] nth number of threads to be used in the given module
330 void put_atoms_in_box_omp(int ePBC, const matrix box, gmx::ArrayRef<gmx::RVec> x, gmx_unused int nth);
332 /*! \brief Put atoms inside triclinic box
334 * This puts ALL atoms in the triclinic unit cell, centered around the
335 * box center as calculated by calc_box_center.
336 * \param[in] ecenter The pbc center type
337 * \param[in] box The simulation box
338 * \param[in,out] x The coordinates of the atoms
340 void put_atoms_in_triclinic_unitcell(int ecenter, const matrix box, gmx::ArrayRef<gmx::RVec> x);
342 /*! \brief Put atoms inside the unitcell
344 * This puts ALL atoms at the closest distance for the center of the box
345 * as calculated by calc_box_center.
346 * When ePBC=-1, the type of pbc is guessed from the box matrix.
347 * \param[in] ePBC The pbc type
348 * \param[in] ecenter The pbc center type
349 * \param[in] box The simulation box
350 * \param[in,out] x The coordinates of the atoms
352 void put_atoms_in_compact_unitcell(int ePBC, int ecenter, const matrix box, gmx::ArrayRef<gmx::RVec> x);
354 /*! \brief Make all molecules whole by shifting positions
356 * \param[in] fplog Log file
357 * \param[in] ePBC The PBC type
358 * \param[in] box The simulation box
359 * \param[in] mtop System topology definition
360 * \param[in,out] x The coordinates of the atoms
362 void do_pbc_first_mtop(FILE* fplog, int ePBC, const matrix box, const gmx_mtop_t* mtop, rvec x[]);
364 /*! \brief Make molecules consisting of multiple charge groups whole by shifting positions
366 * \param[in] ePBC The PBC type
367 * \param[in] box The simulation box
368 * \param[in] mtop System topology definition
369 * \param[in,out] x The coordinates of the atoms
371 void do_pbc_mtop(int ePBC, const matrix box, const gmx_mtop_t* mtop, rvec x[]);