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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * 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.
41 #include "visibility.h"
44 #include "types/commrec.h"
50 #define TRICLINIC(box) (box[YY][XX] != 0 || box[ZZ][XX] != 0 || box[ZZ][YY] != 0)
57 ecenterTRIC
, /* 0.5*(a+b+c) */
58 ecenterRECT
, /* (0.5*a[x],0.5*b[y],0.5*c[z]) */
59 ecenterZERO
, /* (0,0,0) */
60 ecenterDEF
= ecenterTRIC
64 int ePBC2npbcdim(int ePBC
);
65 /* Returns the number of dimensions that use pbc, starting at X */
68 int inputrec2nboundeddim(t_inputrec
*ir
);
69 /* Returns the number of dimensions in which
70 * the coordinates of the particles are bounded, starting at X.
73 void dump_pbc(FILE *fp
, t_pbc
*pbc
);
74 /* Dump the contents of the pbc structure to the file */
77 const char *check_box(int ePBC
, matrix box
);
78 /* Returns NULL if the box is supported by Gromacs.
79 * Otherwise is returns a string with the problem.
80 * When ePBC=-1, the type of pbc is guessed from the box matrix.
84 real
max_cutoff2(int ePBC
, matrix box
);
85 /* Returns the square of the maximum cut-off allowed for the box,
86 * taking into account that the grid neighborsearch code and pbc_dx
87 * only check combinations of single box-vector shifts.
91 int guess_ePBC(matrix box
);
92 /* Guesses the type of periodic boundary conditions using the box */
95 gmx_bool
correct_box(FILE *fplog
, int step
, tensor box
, t_graph
*graph
);
96 /* Checks for un-allowed box angles and corrects the box
97 * and the integer shift vectors in the graph (if graph!=NULL) if necessary.
98 * Returns TRUE when the box was corrected.
102 int ndof_com(t_inputrec
*ir
);
103 /* Returns the number of degrees of freedom of the center of mass */
106 void set_pbc(t_pbc
*pbc
, int ePBC
, matrix box
);
107 /* Initiate the periodic boundary conditions.
108 * pbc_dx will not use pbc and return the normal difference vector
109 * when one or more of the diagonal elements of box are zero.
110 * When ePBC=-1, the type of pbc is guessed from the box matrix.
114 t_pbc
*set_pbc_dd(t_pbc
*pbc
, int ePBC
,
115 gmx_domdec_t
*dd
, gmx_bool bSingleDir
, matrix box
);
116 /* As set_pbc, but additionally sets that correct distances can
117 * be obtained using (combinations of) single box-vector shifts.
118 * Should be used with pbc_dx_aiuc.
119 * If dd!=NULL pbc is not used for directions
120 * with dd->nc[i]==1 with bSingleDir==TRUE or
121 * with dd->nc[i]<=2 with bSingleDir==FALSE.
122 * Returns pbc when pbc operations are required, NULL otherwise.
126 void pbc_dx(const t_pbc
*pbc
, const rvec x1
, const rvec x2
, rvec dx
);
127 /* Calculate the correct distance vector from x2 to x1 and put it in dx.
128 * set_pbc must be called before ever calling this routine.
130 * For triclinic boxes pbc_dx does not necessarily return the shortest
131 * distance vector. If pbc->bLimitDistance=TRUE an atom pair with
132 * distance vector dx with norm2(dx) > pbc->limit_distance2 could
133 * have a shorter distance, but not shorter than sqrt(pbc->limit_distance2).
134 * pbc->limit_distance2 is always larger than max_cutoff2(box).
135 * For the standard rhombic dodecahedron and truncated octahedron
136 * pbc->bLimitDistance=FALSE and thus all distances are correct.
140 int pbc_dx_aiuc(const t_pbc
*pbc
, const rvec x1
, const rvec x2
, rvec dx
);
141 /* Calculate the correct distance vector from x2 to x1 and put it in dx,
142 * This function can only be used when all atoms are in the rectangular
143 * or triclinic unit-cell.
144 * Returns the ishift required to shift x1 at closest distance to x2;
145 * i.e. if 0<=ishift<SHIFTS then x1 - x2 + shift_vec[ishift] = dx
146 * (see calc_shifts below on how to obtain shift_vec)
147 * set_pbc_dd or set_pbc must be called before ever calling this routine.
150 void pbc_dx_d(const t_pbc
*pbc
, const dvec x1
, const dvec x2
, dvec dx
);
151 /* As pbc_dx, but for double precision vectors.
152 * set_pbc must be called before ever calling this routine.
155 gmx_bool
image_rect(ivec xi
, ivec xj
, ivec box_size
,
156 real rlong2
, int *shift
, real
*r2
);
157 /* Calculate the distance between xi and xj for a rectangular box.
158 * When the distance is SMALLER than rlong2 return TRUE, return
159 * the shift code in shift and the distance in r2. When the distance is
160 * >= rlong2 return FALSE;
161 * It is assumed that rlong2 is scaled the same way as the ivecs xi and xj.
164 gmx_bool
image_tri(ivec xi
, ivec xj
, imatrix box
,
165 real rlong2
, int *shift
, real
*r2
);
166 /* Calculate the distance between xi and xj for a triclinic box.
167 * When the distance is SMALLER than rlong2 return TRUE, return
168 * the shift code in shift and the distance in r2. When the distance is
169 * >= rlong2 return FALSE;
170 * It is assumed that rlong2 is scaled the same way as the ivecs xi and xj.
173 gmx_bool
image_cylindric(ivec xi
, ivec xj
, ivec box_size
, real rlong2
,
174 int *shift
, real
*r2
);
175 /* Calculate the distance between xi and xj for a rectangular box
176 * using a cylindric cutoff for long-range only.
177 * When the distance is SMALLER than rlong2 (in X and Y dir.)
178 * return TRUE, return
179 * the shift code in shift and the distance in r2. When the distance is
180 * >= rlong2 return FALSE;
181 * It is assumed that rlong2 is scaled the same way as the ivecs xi and xj.
185 void calc_shifts(matrix box
, rvec shift_vec
[]);
186 /* This routine calculates ths shift vectors necessary to use the
191 void calc_box_center(int ecenter
, matrix box
, rvec box_center
);
192 /* Calculates the center of the box.
193 * See the description for the enum ecenter above.
197 void calc_triclinic_images(matrix box
, rvec img
[]);
198 /* Calculates the NTRICIMG box images */
201 void calc_compact_unitcell_vertices(int ecenter
, matrix box
,
203 /* Calculates the NCUCVERT vertices of a compact unitcell */
206 int *compact_unitcell_edges(void);
207 /* Return an array of unitcell edges of length NCUCEDGE*2,
208 * this is an index in vert[], which is calculated by calc_unitcell_vertices.
209 * The index consists of NCUCEDGE pairs of vertex indices.
210 * The index does not change, so it needs to be retrieved only once.
214 void put_atoms_in_box_omp(int ePBC
, matrix box
, int natoms
, rvec x
[]);
215 /* This wrapper function around put_atoms_in_box() with the ugly manual
216 * workload splitting is needed toavoid silently introducing multithreading
222 void put_atoms_in_box(int ePBC
, matrix box
, int natoms
, rvec x
[]);
223 /* These routines puts ONE or ALL atoms in the box, not caring
224 * about charge groups!
225 * Also works for triclinic cells.
229 void put_atoms_in_triclinic_unitcell(int ecenter
, matrix box
,
230 int natoms
, rvec x
[]);
231 /* This puts ALL atoms in the triclinic unit cell, centered around the
232 * box center as calculated by calc_box_center.
236 const char *put_atoms_in_compact_unitcell(int ePBC
, int ecenter
,
238 int natoms
, rvec x
[]);
239 /* This puts ALL atoms at the closest distance for the center of the box
240 * as calculated by calc_box_center.
241 * Will return NULL is everything went ok and a warning string if not
242 * all atoms could be placed in the unitcell. This can happen for some
243 * triclinic unitcells, see the comment at pbc_dx above.
244 * When ePBC=-1, the type of pbc is guessed from the box matrix.