6efaf9adc1d85d7216b3e88cf71c5e4c3b89f9ec
[alexxy/gromacs.git] / src / gromacs / legacyheaders / pbc.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2012, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37
38 #ifndef _types_pbc_h
39 #define _types_pbc_h
40
41 #include "sysstuff.h"
42 #include "typedefs.h"
43 #include "types/commrec.h"
44
45 #ifdef __cplusplus
46 extern "C" {
47 #endif
48
49 #define TRICLINIC(box) (box[YY][XX] != 0 || box[ZZ][XX] != 0 || box[ZZ][YY] != 0)
50
51 #define NTRICIMG 14
52 #define NCUCVERT 24
53 #define NCUCEDGE 36
54
55 enum {
56     ecenterTRIC, /* 0.5*(a+b+c)                  */
57     ecenterRECT, /* (0.5*a[x],0.5*b[y],0.5*c[z]) */
58     ecenterZERO, /* (0,0,0)                      */
59     ecenterDEF = ecenterTRIC
60 };
61
62 int ePBC2npbcdim(int ePBC);
63 /* Returns the number of dimensions that use pbc, starting at X */
64
65 int inputrec2nboundeddim(t_inputrec *ir);
66 /* Returns the number of dimensions in which
67  * the coordinates of the particles are bounded, starting at X.
68  */
69
70 void dump_pbc(FILE *fp, t_pbc *pbc);
71 /* Dump the contents of the pbc structure to the file */
72
73 const char *check_box(int ePBC, matrix box);
74 /* Returns NULL if the box is supported by Gromacs.
75  * Otherwise is returns a string with the problem.
76  * When ePBC=-1, the type of pbc is guessed from the box matrix.
77  */
78
79 real max_cutoff2(int ePBC, matrix box);
80 /* Returns the square of the maximum cut-off allowed for the box,
81  * taking into account that the grid neighborsearch code and pbc_dx
82  * only check combinations of single box-vector shifts.
83  */
84
85 int guess_ePBC(matrix box);
86 /* Guesses the type of periodic boundary conditions using the box */
87
88 gmx_bool correct_box(FILE *fplog, int step, tensor box, t_graph *graph);
89 /* Checks for un-allowed box angles and corrects the box
90  * and the integer shift vectors in the graph (if graph!=NULL) if necessary.
91  * Returns TRUE when the box was corrected.
92  */
93
94 int ndof_com(t_inputrec *ir);
95 /* Returns the number of degrees of freedom of the center of mass */
96
97 void set_pbc(t_pbc *pbc, int ePBC, matrix box);
98 /* Initiate the periodic boundary conditions.
99  * pbc_dx will not use pbc and return the normal difference vector
100  * when one or more of the diagonal elements of box are zero.
101  * When ePBC=-1, the type of pbc is guessed from the box matrix.
102  */
103
104 t_pbc *set_pbc_dd(t_pbc *pbc, int ePBC,
105                   gmx_domdec_t *dd, gmx_bool bSingleDir, matrix box);
106 /* As set_pbc, but additionally sets that correct distances can
107  * be obtained using (combinations of) single box-vector shifts.
108  * Should be used with pbc_dx_aiuc.
109  * If dd!=NULL pbc is not used for directions
110  * with dd->nc[i]==1 with bSingleDir==TRUE or
111  * with dd->nc[i]<=2 with bSingleDir==FALSE.
112  * Returns pbc when pbc operations are required, NULL otherwise.
113  */
114
115 void pbc_dx(const t_pbc *pbc, const rvec x1, const rvec x2, rvec dx);
116 /* Calculate the correct distance vector from x2 to x1 and put it in dx.
117  * set_pbc must be called before ever calling this routine.
118  *
119  * For triclinic boxes pbc_dx does not necessarily return the shortest
120  * distance vector. If pbc->bLimitDistance=TRUE an atom pair with
121  * distance vector dx with norm2(dx) > pbc->limit_distance2 could
122  * have a shorter distance, but not shorter than sqrt(pbc->limit_distance2).
123  * pbc->limit_distance2 is always larger than max_cutoff2(box).
124  * For the standard rhombic dodecahedron and truncated octahedron
125  * pbc->bLimitDistance=FALSE and thus all distances are correct.
126  */
127
128 int pbc_dx_aiuc(const t_pbc *pbc, const rvec x1, const rvec x2, rvec dx);
129 /* Calculate the correct distance vector from x2 to x1 and put it in dx,
130  * This function can only be used when all atoms are in the rectangular
131  * or triclinic unit-cell.
132  * Returns the ishift required to shift x1 at closest distance to x2;
133  * i.e. if 0<=ishift<SHIFTS then x1 - x2 + shift_vec[ishift] = dx
134  * (see calc_shifts below on how to obtain shift_vec)
135  * set_pbc_dd or set_pbc must be called before ever calling this routine.
136  */
137 void pbc_dx_d(const t_pbc *pbc, const dvec x1, const dvec x2, dvec dx);
138 /* As pbc_dx, but for double precision vectors.
139  * set_pbc must be called before ever calling this routine.
140  */
141
142 gmx_bool image_rect(ivec xi, ivec xj, ivec box_size,
143                     real rlong2, int *shift, real *r2);
144 /* Calculate the distance between xi and xj for a rectangular box.
145  * When the distance is SMALLER than rlong2 return TRUE, return
146  * the shift code in shift and the distance in r2. When the distance is
147  * >= rlong2 return FALSE;
148  * It is assumed that rlong2 is scaled the same way as the ivecs xi and xj.
149  */
150
151 gmx_bool image_tri(ivec xi, ivec xj, imatrix box,
152                    real rlong2, int *shift, real *r2);
153 /* Calculate the distance between xi and xj for a triclinic box.
154  * When the distance is SMALLER than rlong2 return TRUE, return
155  * the shift code in shift and the distance in r2. When the distance is
156  * >= rlong2 return FALSE;
157  * It is assumed that rlong2 is scaled the same way as the ivecs xi and xj.
158  */
159
160 gmx_bool image_cylindric(ivec xi, ivec xj, ivec box_size, real rlong2,
161                          int *shift, real *r2);
162 /* Calculate the distance between xi and xj for a rectangular box
163  * using a cylindric cutoff for long-range only.
164  * When the distance is SMALLER than rlong2 (in X and Y dir.)
165  * return TRUE, return
166  * the shift code in shift and the distance in r2. When the distance is
167  * >= rlong2 return FALSE;
168  * It is assumed that rlong2 is scaled the same way as the ivecs xi and xj.
169  */
170
171 void calc_shifts(matrix box, rvec shift_vec[]);
172 /* This routine calculates ths shift vectors necessary to use the
173  * ns routine.
174  */
175
176 void calc_box_center(int ecenter, matrix box, rvec box_center);
177 /* Calculates the center of the box.
178  * See the description for the enum ecenter above.
179  */
180
181 void calc_triclinic_images(matrix box, rvec img[]);
182 /* Calculates the NTRICIMG box images */
183
184 void calc_compact_unitcell_vertices(int ecenter, matrix box,
185                                     rvec vert[]);
186 /* Calculates the NCUCVERT vertices of a compact unitcell */
187
188 int *compact_unitcell_edges(void);
189 /* Return an array of unitcell edges of length NCUCEDGE*2,
190  * this is an index in vert[], which is calculated by calc_unitcell_vertices.
191  * The index consists of NCUCEDGE pairs of vertex indices.
192  * The index does not change, so it needs to be retrieved only once.
193  */
194
195 void put_atoms_in_box_omp(int ePBC, matrix box, int natoms, rvec x[]);
196 /* This wrapper function around put_atoms_in_box() with the ugly manual
197  * workload splitting is needed toavoid silently introducing multithreading
198  * in tools.
199  * */
200
201
202 void put_atoms_in_box(int ePBC, matrix box, int natoms, rvec x[]);
203 /* These routines puts ONE or ALL atoms in the box, not caring
204  * about charge groups!
205  * Also works for triclinic cells.
206  */
207
208 void put_atoms_in_triclinic_unitcell(int ecenter, matrix box,
209                                      int natoms, rvec x[]);
210 /* This puts ALL atoms in the triclinic unit cell, centered around the
211  * box center as calculated by calc_box_center.
212  */
213
214 const char *put_atoms_in_compact_unitcell(int ePBC, int ecenter,
215                                           matrix box,
216                                           int natoms, rvec x[]);
217 /* This puts ALL atoms at the closest distance for the center of the box
218  * as calculated by calc_box_center.
219  * Will return NULL is everything went ok and a warning string if not
220  * all atoms could be placed in the unitcell. This can happen for some
221  * triclinic unitcells, see the comment at pbc_dx above.
222  * When ePBC=-1, the type of pbc is guessed from the box matrix.
223  */
224
225 #ifdef __cplusplus
226 }
227 #endif
228
229 #endif  /* _pbc_h */