added Verlet scheme and NxN non-bonded functionality
[gromacs.git] / include / pbc.h
blobe3c49ca249a2ad0f6e9e9553293ca872887fe51b
1 /*
2 *
3 * This source code is part of
4 *
5 * G R O M A C S
6 *
7 * GROningen MAchine for Chemical Simulations
8 *
9 * VERSION 3.2.0
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
32 * And Hey:
33 * Gromacs Runs On Most of All Computer Systems
36 #ifndef _types_pbc_h
37 #define _types_pbc_h
39 #include "sysstuff.h"
40 #include "typedefs.h"
41 #include "types/commrec.h"
43 #ifdef __cplusplus
44 extern "C" {
45 #endif
47 #define TRICLINIC(box) (box[YY][XX]!=0 || box[ZZ][XX]!=0 || box[ZZ][YY]!=0)
49 #define NTRICIMG 14
50 #define NCUCVERT 24
51 #define NCUCEDGE 36
53 enum {
54 ecenterTRIC, /* 0.5*(a+b+c) */
55 ecenterRECT, /* (0.5*a[x],0.5*b[y],0.5*c[z]) */
56 ecenterZERO, /* (0,0,0) */
57 ecenterDEF = ecenterTRIC
60 int ePBC2npbcdim(int ePBC);
61 /* Returns the number of dimensions that use pbc, starting at X */
63 int inputrec2nboundeddim(t_inputrec *ir);
64 /* Returns the number of dimensions in which
65 * the coordinates of the particles are bounded, starting at X.
68 void dump_pbc(FILE *fp,t_pbc *pbc);
69 /* Dump the contents of the pbc structure to the file */
71 const char *check_box(int ePBC,matrix box);
72 /* Returns NULL if the box is supported by Gromacs.
73 * Otherwise is returns a string with the problem.
74 * When ePBC=-1, the type of pbc is guessed from the box matrix.
77 real max_cutoff2(int ePBC,matrix box);
78 /* Returns the square of the maximum cut-off allowed for the box,
79 * taking into account that the grid neighborsearch code and pbc_dx
80 * only check combinations of single box-vector shifts.
83 int guess_ePBC(matrix box);
84 /* Guesses the type of periodic boundary conditions using the box */
86 gmx_bool correct_box(FILE *fplog,int step,tensor box,t_graph *graph);
87 /* Checks for un-allowed box angles and corrects the box
88 * and the integer shift vectors in the graph (if graph!=NULL) if necessary.
89 * Returns TRUE when the box was corrected.
92 int ndof_com(t_inputrec *ir);
93 /* Returns the number of degrees of freedom of the center of mass */
95 void set_pbc(t_pbc *pbc,int ePBC,matrix box);
96 /* Initiate the periodic boundary conditions.
97 * pbc_dx will not use pbc and return the normal difference vector
98 * when one or more of the diagonal elements of box are zero.
99 * When ePBC=-1, the type of pbc is guessed from the box matrix.
102 t_pbc *set_pbc_dd(t_pbc *pbc,int ePBC,
103 gmx_domdec_t *dd,gmx_bool bSingleDir,matrix box);
104 /* As set_pbc, but additionally sets that correct distances can
105 * be obtained using (combinations of) single box-vector shifts.
106 * Should be used with pbc_dx_aiuc.
107 * If dd!=NULL pbc is not used for directions
108 * with dd->nc[i]==1 with bSingleDir==TRUE or
109 * with dd->nc[i]<=2 with bSingleDir==FALSE.
110 * Returns pbc when pbc operations are required, NULL otherwise.
113 void pbc_dx(const t_pbc *pbc,const rvec x1, const rvec x2, rvec dx);
114 /* Calculate the correct distance vector from x2 to x1 and put it in dx.
115 * set_pbc must be called before ever calling this routine.
117 * For triclinic boxes pbc_dx does not necessarily return the shortest
118 * distance vector. If pbc->bLimitDistance=TRUE an atom pair with
119 * distance vector dx with norm2(dx) > pbc->limit_distance2 could
120 * have a shorter distance, but not shorter than sqrt(pbc->limit_distance2).
121 * pbc->limit_distance2 is always larger than max_cutoff2(box).
122 * For the standard rhombic dodecahedron and truncated octahedron
123 * pbc->bLimitDistance=FALSE and thus all distances are correct.
126 int pbc_dx_aiuc(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 * This function can only be used when all atoms are in the rectangular
129 * or triclinic unit-cell.
130 * Returns the ishift required to shift x1 at closest distance to x2;
131 * i.e. if 0<=ishift<SHIFTS then x1 - x2 + shift_vec[ishift] = dx
132 * (see calc_shifts below on how to obtain shift_vec)
133 * set_pbc_dd or set_pbc must be called before ever calling this routine.
135 void pbc_dx_d(const t_pbc *pbc,const dvec x1, const dvec x2, dvec dx);
136 /* As pbc_dx, but for double precision vectors.
137 * set_pbc must be called before ever calling this routine.
140 gmx_bool image_rect(ivec xi,ivec xj,ivec box_size,
141 real rlong2,int *shift,real *r2);
142 /* Calculate the distance between xi and xj for a rectangular box.
143 * When the distance is SMALLER than rlong2 return TRUE, return
144 * the shift code in shift and the distance in r2. When the distance is
145 * >= rlong2 return FALSE;
146 * It is assumed that rlong2 is scaled the same way as the ivecs xi and xj.
149 gmx_bool image_tri(ivec xi,ivec xj,imatrix box,
150 real rlong2,int *shift,real *r2);
151 /* Calculate the distance between xi and xj for a triclinic box.
152 * When the distance is SMALLER than rlong2 return TRUE, return
153 * the shift code in shift and the distance in r2. When the distance is
154 * >= rlong2 return FALSE;
155 * It is assumed that rlong2 is scaled the same way as the ivecs xi and xj.
158 gmx_bool image_cylindric(ivec xi,ivec xj,ivec box_size,real rlong2,
159 int *shift,real *r2);
160 /* Calculate the distance between xi and xj for a rectangular box
161 * using a cylindric cutoff for long-range only.
162 * When the distance is SMALLER than rlong2 (in X and Y dir.)
163 * return TRUE, return
164 * the shift code in shift and the distance in r2. When the distance is
165 * >= rlong2 return FALSE;
166 * It is assumed that rlong2 is scaled the same way as the ivecs xi and xj.
169 void calc_shifts(matrix box,rvec shift_vec[]);
170 /* This routine calculates ths shift vectors necessary to use the
171 * ns routine.
174 void calc_box_center(int ecenter,matrix box,rvec box_center);
175 /* Calculates the center of the box.
176 * See the description for the enum ecenter above.
179 void calc_triclinic_images(matrix box,rvec img[]);
180 /* Calculates the NTRICIMG box images */
182 void calc_compact_unitcell_vertices(int ecenter,matrix box,
183 rvec vert[]);
184 /* Calculates the NCUCVERT vertices of a compact unitcell */
186 int *compact_unitcell_edges(void);
187 /* Return an array of unitcell edges of length NCUCEDGE*2,
188 * this is an index in vert[], which is calculated by calc_unitcell_vertices.
189 * The index consists of NCUCEDGE pairs of vertex indices.
190 * The index does not change, so it needs to be retrieved only once.
193 void put_atoms_in_box_omp(int ePBC,matrix box,int natoms,rvec x[]);
194 /* This wrapper function around put_atoms_in_box() with the ugly manual
195 * workload splitting is needed toavoid silently introducing multithreading
196 * in tools.
197 * */
200 void put_atoms_in_box(int ePBC, matrix box,int natoms,rvec x[]);
201 /* These routines puts ONE or ALL atoms in the box, not caring
202 * about charge groups!
203 * Also works for triclinic cells.
206 void put_atoms_in_triclinic_unitcell(int ecenter,matrix box,
207 int natoms,rvec x[]);
208 /* This puts ALL atoms in the triclinic unit cell, centered around the
209 * box center as calculated by calc_box_center.
212 const char *put_atoms_in_compact_unitcell(int ePBC,int ecenter,
213 matrix box,
214 int natoms,rvec x[]);
215 /* This puts ALL atoms at the closest distance for the center of the box
216 * as calculated by calc_box_center.
217 * Will return NULL is everything went ok and a warning string if not
218 * all atoms could be placed in the unitcell. This can happen for some
219 * triclinic unitcells, see the comment at pbc_dx above.
220 * When ePBC=-1, the type of pbc is guessed from the box matrix.
223 #ifdef __cplusplus
225 #endif
227 #endif /* _pbc_h */