3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
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
33 * Gromacs Runs On Most of All Computer Systems
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
63 extern int ePBC2npbcdim(int ePBC
);
64 /* Returns the number of dimensions that use pbc, starting at X */
66 extern int inputrec2nboundeddim(t_inputrec
*ir
);
67 /* Returns the number of dimensions in which
68 * the coordinates of the particles are bounded, starting at X.
71 extern void dump_pbc(FILE *fp
,t_pbc
*pbc
);
72 /* Dump the contents of the pbc structure to the file */
74 extern const char *check_box(int ePBC
,matrix box
);
75 /* Returns NULL if the box is supported by Gromacs.
76 * Otherwise is returns a string with the problem.
77 * When ePBC=-1, the type of pbc is guessed from the box matrix.
80 extern real
max_cutoff2(int ePBC
,matrix box
);
81 /* Returns the square of the maximum cut-off allowed for the box,
82 * taking into account that the grid neighborsearch code and pbc_dx
83 * only check combinations of single box-vector shifts.
86 int guess_ePBC(matrix box
);
87 /* Guesses the type of periodic boundary conditions using the box */
89 extern bool correct_box(FILE *fplog
,int step
,tensor box
,t_graph
*graph
);
90 /* Checks for un-allowed box angles and corrects the box
91 * and the integer shift vectors in the graph (if graph!=NULL) if necessary.
92 * Returns TRUE when the box was corrected.
95 extern int ndof_com(t_inputrec
*ir
);
96 /* Returns the number of degrees of freedom of the center of mass */
98 extern void set_pbc(t_pbc
*pbc
,int ePBC
,matrix box
);
99 /* Initiate the periodic boundary conditions.
100 * pbc_dx will not use pbc and return the normal difference vector
101 * when one or more of the diagonal elements of box are zero.
102 * When ePBC=-1, the type of pbc is guessed from the box matrix.
105 extern t_pbc
*set_pbc_dd(t_pbc
*pbc
,int ePBC
,
106 gmx_domdec_t
*dd
,bool bSingleDir
,matrix box
);
107 /* As set_pbc, but additionally sets that correct distances can
108 * be obtained using (combinations of) single box-vector shifts.
109 * Should be used with pbc_dx_aiuc.
110 * If dd!=NULL pbc is not used for directions
111 * with dd->nc[i]==1 with bSingleDir==TRUE or
112 * with dd->nc[i]<=2 with bSingleDir==FALSE.
113 * Returns pbc when pbc operations are required, NULL otherwise.
116 extern void pbc_dx(const t_pbc
*pbc
,const rvec x1
, const rvec x2
, rvec dx
);
117 /* Calculate the correct distance vector from x2 to x1 and put it in dx.
118 * set_pbc must be called before ever calling this routine.
120 * For triclinic boxes pbc_dx does not necessarily return the shortest
121 * distance vector. If pbc->bLimitDistance=TRUE an atom pair with
122 * distance vector dx with norm2(dx) > pbc->limit_distance2 could
123 * have a shorter distance, but not shorter than sqrt(pbc->limit_distance2).
124 * pbc->limit_distance2 is always larger than max_cutoff2(box).
125 * For the standard rhombic dodecahedron and truncated octahedron
126 * pbc->bLimitDistance=FALSE and thus all distances are correct.
129 extern int pbc_dx_aiuc(const t_pbc
*pbc
,const rvec x1
,const rvec x2
,rvec dx
);
130 /* Calculate the correct distance vector from x2 to x1 and put it in dx,
131 * This function can only be used when all atoms are in the rectangular
132 * or triclinic unit-cell.
133 * Returns the ishift required to shift x1 at closest distance to x2;
134 * i.e. if 0<=ishift<SHIFTS then x1 - x2 + shift_vec[ishift] = dx
135 * (see calc_shifts below on how to obtain shift_vec)
136 * set_pbc_dd or set_pbc must be called before ever calling this routine.
138 extern void pbc_dx_d(const t_pbc
*pbc
,const dvec x1
, const dvec x2
, dvec dx
);
139 /* As pbc_dx, but for double precision vectors.
140 * set_pbc must be called before ever calling this routine.
143 extern bool image_rect(ivec xi
,ivec xj
,ivec box_size
,
144 real rlong2
,int *shift
,real
*r2
);
145 /* Calculate the distance between xi and xj for a rectangular box.
146 * When the distance is SMALLER than rlong2 return TRUE, return
147 * the shift code in shift and the distance in r2. When the distance is
148 * >= rlong2 return FALSE;
149 * It is assumed that rlong2 is scaled the same way as the ivecs xi and xj.
152 extern bool image_tri(ivec xi
,ivec xj
,imatrix box
,
153 real rlong2
,int *shift
,real
*r2
);
154 /* Calculate the distance between xi and xj for a triclinic box.
155 * When the distance is SMALLER than rlong2 return TRUE, return
156 * the shift code in shift and the distance in r2. When the distance is
157 * >= rlong2 return FALSE;
158 * It is assumed that rlong2 is scaled the same way as the ivecs xi and xj.
161 extern bool image_cylindric(ivec xi
,ivec xj
,ivec box_size
,real rlong2
,
162 int *shift
,real
*r2
);
163 /* Calculate the distance between xi and xj for a rectangular box
164 * using a cylindric cutoff for long-range only.
165 * When the distance is SMALLER than rlong2 (in X and Y dir.)
166 * return TRUE, return
167 * the shift code in shift and the distance in r2. When the distance is
168 * >= rlong2 return FALSE;
169 * It is assumed that rlong2 is scaled the same way as the ivecs xi and xj.
172 extern void calc_shifts(matrix box
,rvec shift_vec
[]);
173 /* This routine calculates ths shift vectors necessary to use the
177 extern void calc_cgcm(FILE *log
,int cg0
,int cg1
,t_block
*cgs
,
178 rvec pos
[],rvec cg_cm
[]);
179 /* Routine to compute centers of geometry of charge groups. No periodicity
183 extern void put_charge_groups_in_box (FILE *log
,int cg0
,int cg1
,
184 int ePBC
,matrix box
,t_block
*cgs
,
188 /* This routine puts charge groups in the periodic box, keeping them
192 extern void calc_box_center(int ecenter
,matrix box
,rvec box_center
);
193 /* Calculates the center of the box.
194 * See the description for the enum ecenter above.
197 extern void calc_triclinic_images(matrix box
,rvec img
[]);
198 /* Calculates the NTRICIMG box images */
200 extern void calc_compact_unitcell_vertices(int ecenter
,matrix box
,
202 /* Calculates the NCUCVERT vertices of a compact unitcell */
204 extern int *compact_unitcell_edges(void);
205 /* Return an array of unitcell edges of length NCUCEDGE*2,
206 * this is an index in vert[], which is calculated by calc_unitcell_vertices.
207 * The index consists of NCUCEDGE pairs of vertex indices.
208 * The index does not change, so it needs to be retrieved only once.
210 extern void put_atom_in_box(matrix box
,rvec x
);
212 extern void put_atoms_in_box(matrix box
,int natoms
,rvec x
[]);
213 /* These routines puts ONE or ALL atoms in the box, not caring
214 * about charge groups!
215 * Also works for triclinic cells.
218 extern void put_atoms_in_triclinic_unitcell(int ecenter
,matrix box
,
219 int natoms
,rvec x
[]);
220 /* This puts ALL atoms in the triclinic unit cell, centered around the
221 * box center as calculated by calc_box_center.
224 extern const char *put_atoms_in_compact_unitcell(int ePBC
,int ecenter
,
226 int natoms
,rvec x
[]);
227 /* This puts ALL atoms at the closest distance for the center of the box
228 * as calculated by calc_box_center.
229 * Will return NULL is everything went ok and a warning string if not
230 * all atoms could be placed in the unitcell. This can happen for some
231 * triclinic unitcells, see the comment at pbc_dx above.
232 * When ePBC=-1, the type of pbc is guessed from the box matrix.