Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / pbcutil / pbc.h
blobdc70f22994e26a4cfc44d2f1d861dbb794f0d959
1 /*
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,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 #ifndef GMX_PBCUTIL_PBC_H
38 #define GMX_PBCUTIL_PBC_H
40 #include <stdio.h>
42 #include "gromacs/math/vectypes.h"
43 #include "gromacs/utility/arrayref.h"
44 #include "gromacs/utility/basedefinitions.h"
45 #include "gromacs/utility/real.h"
47 struct gmx_domdec_t;
48 struct gmx_mtop_t;
50 enum {
51 epbcXYZ, epbcNONE, epbcXY, epbcSCREW, epbcNR
54 //! Strings corresponding to epbc enum values.
55 extern const char *epbc_names[epbcNR+1];
57 /* Maximum number of combinations of single triclinic box vectors
58 * required to shift atoms that are within a brick of the size of
59 * the diagonal of the box to within the maximum cut-off distance.
61 #define MAX_NTRICVEC 12
63 /*! \brief Structure containing info on periodic boundary conditions */
64 typedef struct t_pbc {
65 //! The PBC type
66 int ePBC;
67 //! Number of dimensions in which PBC is exerted
68 int ndim_ePBC;
69 /*! \brief Determines how to compute distance vectors.
71 * Indicator of how to compute distance vectors, depending
72 * on PBC type (depends on ePBC and dimensions with(out) DD)
73 * and the box angles.
75 int ePBCDX;
76 /*! \brief Used for selecting which dimensions to use in PBC.
78 * In case of 1-D PBC this indicates which dimension is used,
79 * in case of 2-D PBC this indicates the opposite
81 int dim;
82 //! The simulation box
83 matrix box;
84 //! The lengths of the diagonal of the full box
85 rvec fbox_diag;
86 //! Halve of the above
87 rvec hbox_diag;
88 //! Negative of the above
89 rvec mhbox_diag;
90 //! Maximum allowed cutoff squared for the box and PBC used
91 real max_cutoff2;
92 /*! \brief Number of triclinic shift vectors.
94 * Number of triclinic shift vectors depends on the skewedness
95 * of the box, that is mostly on the angles. For triclinic boxes
96 * we first take the closest image along each Cartesian dimension
97 * independently. When the resulting distance^2 is larger than
98 * max_cutoff2, up to ntric_vec triclinic shift vectors need to
99 * be tried. Because of the restrictions imposed on the unit-cell
100 * by GROMACS, ntric_vec <= MAX_NTRICVEC = 12.
102 int ntric_vec;
103 //! The triclinic shift vectors in grid cells. Internal use only.
104 ivec tric_shift[MAX_NTRICVEC];
105 //! The triclinic shift vectors in length units
106 rvec tric_vec[MAX_NTRICVEC];
107 } t_pbc;
109 #define TRICLINIC(box) ((box)[YY][XX] != 0 || (box)[ZZ][XX] != 0 || (box)[ZZ][YY] != 0)
111 #define NTRICIMG 14
112 #define NCUCVERT 24
113 #define NCUCEDGE 36
115 enum {
116 ecenterTRIC, /* 0.5*(a+b+c) */
117 ecenterRECT, /* (0.5*a[x],0.5*b[y],0.5*c[z]) */
118 ecenterZERO, /* (0,0,0) */
119 ecenterDEF = ecenterTRIC
122 struct t_graph;
124 /*! \brief Returns the number of dimensions that use pbc
126 * \param[in] ePBC The periodic boundary condition type
127 * \return the number of dimensions that use pbc, starting at X
129 int ePBC2npbcdim(int ePBC);
131 /*! \brief Dump the contents of the pbc structure to the file
133 * \param[in] fp The file pointer to write to
134 * \param[in] pbc The periodic boundary condition information structure
136 void dump_pbc(FILE *fp, t_pbc *pbc);
138 /*! \brief Check the box for consistency
140 * \param[in] ePBC The pbc identifier
141 * \param[in] box The box matrix
142 * \return NULL if the box is supported by Gromacs.
143 * Otherwise returns a string with the problem.
144 * When ePBC=-1, the type of pbc is guessed from the box matrix.
146 const char *check_box(int ePBC, const matrix box);
148 /*! \brief Creates box matrix from edge lengths and angles.
150 * \param[in,out] box The box matrix
151 * \param[in] vec The edge lengths
152 * \param[in] angleInDegrees The angles
154 void matrix_convert(matrix box, const rvec vec, const rvec angleInDegrees);
156 /*! \brief Compute the maximum cutoff for the box
158 * Returns the square of the maximum cut-off allowed for the box,
159 * taking into account that the grid neighborsearch code and pbc_dx
160 * only check combinations of single box-vector shifts.
161 * \param[in] ePBC The pbc identifier
162 * \param[in] box The box matrix
163 * \return the maximum cut-off.
165 real max_cutoff2(int ePBC, const matrix box);
167 /*! \brief Guess PBC typr
169 * Guesses the type of periodic boundary conditions using the box
170 * \param[in] box The box matrix
171 * \return The pbc identifier
173 int guess_ePBC(const matrix box);
175 /*! \brief Corrects the box if necessary
177 * Checks for un-allowed box angles and corrects the box
178 * and the integer shift vectors in the graph (if graph!=NULL) if necessary.
179 * \param[in] fplog File for debug output
180 * \param[in] step The MD step number
181 * \param[in] box The simulation cell
182 * \param[in] graph Information about molecular connectivity
183 * \return TRUE when the box was corrected.
185 gmx_bool correct_box(FILE *fplog, int step, tensor box, struct t_graph *graph);
187 /*! \brief Initiate the periodic boundary condition algorithms.
189 * pbc_dx will not use pbc and return the normal difference vector
190 * when one or more of the diagonal elements of box are zero.
191 * When ePBC=-1, the type of pbc is guessed from the box matrix.
192 * \param[in,out] pbc The pbc information structure
193 * \param[in] ePBC The PBC identifier
194 * \param[in] box The box tensor
196 void set_pbc(t_pbc *pbc, int ePBC, const matrix box);
198 /*! \brief Initiate the periodic boundary condition algorithms.
200 * As set_pbc, but additionally sets that correct distances can
201 * be obtained using (combinations of) single box-vector shifts.
202 * Should be used with pbc_dx_aiuc.
203 * If domdecCells!=NULL pbc is not used for directions
204 * with dd->nc[i]==1 with bSingleDir==TRUE or
205 * with dd->nc[i]<=2 with bSingleDir==FALSE.
206 * Note that when no PBC is required only pbc->ePBC is set,
207 * the rest of the struct will be invalid.
208 * \param[in,out] pbc The pbc information structure
209 * \param[in] ePBC The PBC identifier
210 * \param[in] domdecCells 3D integer vector describing the number of DD cells
211 * or nullptr if not using DD.
212 * \param[in] bSingleDir TRUE if DD communicates only in one direction along dimensions
213 * \param[in] box The box tensor
214 * \return the pbc structure when pbc operations are required, NULL otherwise.
216 t_pbc *set_pbc_dd(t_pbc *pbc, int ePBC,
217 const ivec domdecCells, gmx_bool bSingleDir,
218 const matrix box);
220 /*! \brief Compute distance with PBC
222 * Calculate the correct distance vector from x2 to x1 and put it in dx.
223 * set_pbc must be called before ever calling this routine.
225 * Note that for triclinic boxes that do not obey the GROMACS unit-cell
226 * restrictions, pbc_dx and pbc_dx_aiuc will not correct for PBC.
227 * \param[in,out] pbc The pbc information structure
228 * \param[in] x1 Coordinates for particle 1
229 * \param[in] x2 Coordinates for particle 2
230 * \param[out] dx Distance vector
232 void pbc_dx(const t_pbc *pbc, const rvec x1, const rvec x2, rvec dx);
234 /*! \brief Compute distance vector for simple PBC types
236 * Calculate the correct distance vector from x2 to x1 and put it in dx,
237 * This function can only be used when all atoms are in the rectangular
238 * or triclinic unit-cell.
239 * set_pbc_dd or set_pbc must be called before ever calling this routine.
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
244 * \return the ishift required to shift x1 at closest distance to x2;
245 * i.e. if 0<=ishift<SHIFTS then x1 - x2 + shift_vec[ishift] = dx
246 * (see calc_shifts below on how to obtain shift_vec)
248 int pbc_dx_aiuc(const t_pbc *pbc, const rvec x1, const rvec x2, rvec dx);
250 /*! \brief Compute distance with PBC
252 * As pbc_dx, but for double precision vectors.
253 * set_pbc must be called before ever calling this routine.
254 * \param[in,out] pbc The pbc information structure
255 * \param[in] x1 Coordinates for particle 1
256 * \param[in] x2 Coordinates for particle 2
257 * \param[out] dx Distance vector
259 void pbc_dx_d(const t_pbc *pbc, const dvec x1, const dvec x2, dvec dx);
261 /*! \brief Computes shift vectors
263 * This routine calculates ths shift vectors necessary to use the
264 * neighbor searching routine.
265 * \param[in] box The simulation box
266 * \param[out] shift_vec The shifting vectors
268 void calc_shifts(const matrix box, rvec shift_vec[]);
270 /*! \brief Calculates the center of the box.
272 * See the description for the enum ecenter above.
273 * \param[in] ecenter Description of center type
274 * \param[in] box The simulation box
275 * \param[out] box_center The center of the box
277 void calc_box_center(int ecenter, const matrix box, rvec box_center);
279 /*! \brief Calculates the NTRICIMG box images
281 * \param[in] box The simulation box
282 * \param[out] img The triclinic box images
284 void calc_triclinic_images(const matrix box, rvec img[]);
286 /*! \brief Calculates the NCUCVERT vertices of a compact unitcell
288 * \param[in] ecenter The center type
289 * \param[in] box The simulation box
290 * \param[out] vert The vertices
292 void calc_compact_unitcell_vertices(int ecenter, const matrix box,
293 rvec vert[]);
295 /*! \brief Compute unitcell edges
297 * \return an array of unitcell edges of length NCUCEDGE*2,
298 * this is an index in vert[], which is calculated by calc_unitcell_vertices.
299 * The index consists of NCUCEDGE pairs of vertex indices.
300 * The index does not change, so it needs to be retrieved only once.
302 int *compact_unitcell_edges();
304 /*! \brief Put atoms inside the simulations box
306 * These routines puts ONE or ALL atoms in the box, not caring
307 * about charge groups!
308 * Also works for triclinic cells.
309 * \param[in] ePBC The pbc type
310 * \param[in] box The simulation box
311 * \param[in,out] x The coordinates of the atoms
313 void put_atoms_in_box(int ePBC, const matrix box, gmx::ArrayRef<gmx::RVec> x);
315 /*! \brief Parallellizes put_atoms_in_box()
317 * This wrapper function around put_atoms_in_box() with the ugly manual
318 * workload splitting is needed to avoid silently introducing multithreading
319 * in tools.
320 * \param[in] ePBC The pbc type
321 * \param[in] box The simulation box
322 * \param[in,out] x The coordinates of the atoms
323 * \param[in] nth number of threads to be used in the given module
325 void put_atoms_in_box_omp(int ePBC, const matrix box, gmx::ArrayRef<gmx::RVec> x, gmx_unused int nth);
327 /*! \brief Put atoms inside triclinic box
329 * This puts ALL atoms in the triclinic unit cell, centered around the
330 * box center as calculated by calc_box_center.
331 * \param[in] ecenter The pbc center type
332 * \param[in] box The simulation box
333 * \param[in,out] x The coordinates of the atoms
335 void put_atoms_in_triclinic_unitcell(int ecenter, const matrix box,
336 gmx::ArrayRef<gmx::RVec> x);
338 /*! \brief Put atoms inside the unitcell
340 * This puts ALL atoms at the closest distance for the center of the box
341 * as calculated by calc_box_center.
342 * When ePBC=-1, the type of pbc is guessed from the box matrix.
343 * \param[in] ePBC The pbc type
344 * \param[in] ecenter The pbc center type
345 * \param[in] box The simulation box
346 * \param[in,out] x The coordinates of the atoms
348 void put_atoms_in_compact_unitcell(int ePBC, int ecenter,
349 const matrix box,
350 gmx::ArrayRef<gmx::RVec> x);
352 /*! \brief Make all molecules whole by shifting positions
354 * \param[in] fplog Log file
355 * \param[in] ePBC The PBC type
356 * \param[in] box The simulation box
357 * \param[in] mtop System topology definition
358 * \param[in,out] x The coordinates of the atoms
360 void do_pbc_first_mtop(FILE *fplog, int ePBC, const matrix box,
361 const gmx_mtop_t *mtop, rvec x[]);
363 /*! \brief Make molecules consisting of multiple charge groups whole by shifting positions
365 * \param[in] ePBC The PBC type
366 * \param[in] box The simulation box
367 * \param[in] mtop System topology definition
368 * \param[in,out] x The coordinates of the atoms
370 void do_pbc_mtop(int ePBC, const matrix box,
371 const gmx_mtop_t *mtop, rvec x[]);
373 #endif