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