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 * Copyright (c) 2013,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.
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.
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.
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.
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.
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.
42 #include "gromacs/commandline/pargs.h"
43 #include "gromacs/fileio/trxio.h"
44 #include "gromacs/fileio/xvgr.h"
45 #include "gromacs/gmxana/gmx_ana.h"
46 #include "gromacs/gmxana/princ.h"
47 #include "gromacs/legacyheaders/macros.h"
48 #include "gromacs/legacyheaders/typedefs.h"
49 #include "gromacs/legacyheaders/viewit.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/pbcutil/rmpbc.h"
52 #include "gromacs/topology/index.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/futil.h"
55 #include "gromacs/utility/smalloc.h"
57 /****************************************************************************/
58 /* This program calculates the ordering of water molecules across a box, as */
59 /* function of the z-coordinate. This implies averaging over slices and over*/
60 /* time. Output is the average cos of the angle of the dipole with the */
61 /* normal, per slice. */
62 /* In addition, it calculates the average dipole moment itself in three */
64 /****************************************************************************/
66 void calc_h2order(const char *fn, atom_id index[], int ngx, rvec **slDipole,
67 real **slOrder, real *slWidth, int *nslices,
68 t_topology *top, int ePBC,
69 int axis, gmx_bool bMicel, atom_id micel[], int nmic,
70 const output_env_t oenv)
72 rvec *x0, /* coordinates with pbc */
73 dipole, /* dipole moment due to one molecules */
75 com; /* center of mass of micel, with bMicel */
76 rvec *dip; /* sum of dipoles, unnormalized */
77 matrix box; /* box (3x3) */
79 real t, /* time from trajectory */
80 *sum, /* sum of all cosines of dipoles, per slice */
81 *frame; /* order over one frame */
82 int natoms, /* nr. atoms in trj */
84 slice = 0, /* current slice number */
85 *count; /* nr. of atoms in one slice */
86 gmx_rmpbc_t gpbc = NULL;
88 if ((natoms = read_first_x(oenv, &status, fn, &t, &x0, box)) == 0)
90 gmx_fatal(FARGS, "Could not read coordinates from statusfile\n");
95 *nslices = (int)(box[axis][axis] * 10); /* default value */
102 normal[0] = 1; normal[1] = 0; normal[2] = 0;
105 normal[0] = 0; normal[1] = 1; normal[2] = 0;
108 normal[0] = 0; normal[1] = 0; normal[2] = 1;
111 gmx_fatal(FARGS, "No valid value for -axis-. Exiting.\n");
112 /* make compiler happy */
113 normal[0] = 1; normal[1] = 0; normal[2] = 0;
117 snew(count, *nslices);
120 snew(frame, *nslices);
122 *slWidth = box[axis][axis]/(*nslices);
123 fprintf(stderr, "Box divided in %d slices. Initial width of slice: %f\n",
128 gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms);
129 /*********** Start processing trajectory ***********/
132 *slWidth = box[axis][axis]/(*nslices);
135 gmx_rmpbc(gpbc, natoms, box, x0);
139 calc_xcm(x0, nmic, micel, top->atoms.atom, com, FALSE);
142 for (i = 0; i < ngx/3; i++)
144 /* put all waters in box */
145 for (j = 0; j < DIM; j++)
147 if (x0[index[3*i]][j] < 0)
149 x0[index[3*i]][j] += box[j][j];
150 x0[index[3*i+1]][j] += box[j][j];
151 x0[index[3*i+2]][j] += box[j][j];
153 if (x0[index[3*i]][j] > box[j][j])
155 x0[index[3*i]][j] -= box[j][j];
156 x0[index[3*i+1]][j] -= box[j][j];
157 x0[index[3*i+2]][j] -= box[j][j];
161 for (j = 0; j < DIM; j++)
164 x0[index[3*i]][j] * top->atoms.atom[index[3*i]].q +
165 x0[index[3*i+1]][j] * top->atoms.atom[index[3*i+1]].q +
166 x0[index[3*i+2]][j] * top->atoms.atom[index[3*i+2]].q;
169 /* now we have a dipole vector. Might as well safe it. Then the
170 rest depends on whether we're dealing with
171 a flat or a spherical interface.
175 { /* this is for spherical interfaces */
176 rvec_sub(com, x0[index[3*i]], normal); /* vector from Oxygen to COM */
177 slice = norm(normal)/(*slWidth); /* spherical slice */
179 sum[slice] += iprod(dipole, normal) / (norm(dipole) * norm(normal));
180 frame[slice] += iprod(dipole, normal) / (norm(dipole) * norm(normal));
186 /* this is for flat interfaces */
188 /* determine which slice atom is in */
189 slice = (x0[index[3*i]][axis] / (*slWidth));
190 if (slice < 0 || slice >= *nslices)
192 fprintf(stderr, "Coordinate: %f ", x0[index[3*i]][axis]);
193 fprintf(stderr, "HELP PANIC! slice = %d, OUT OF RANGE!\n", slice);
197 rvec_add(dipole, dip[slice], dip[slice]);
198 /* Add dipole to total. mag[slice] is total dipole in axis direction */
199 sum[slice] += iprod(dipole, normal)/norm(dipole);
200 frame[slice] += iprod(dipole, normal)/norm(dipole);
201 /* increase count for that slice */
208 while (read_next_x(oenv, status, &t, x0, box));
209 /*********** done with status file **********/
211 fprintf(stderr, "\nRead trajectory. Printing parameters to file\n");
212 gmx_rmpbc_done(gpbc);
214 for (i = 0; i < *nslices; i++) /* average over frames */
216 fprintf(stderr, "%d waters in slice %d\n", count[i], i);
217 if (count[i] > 0) /* divide by number of molecules in each slice */
219 sum[i] = sum[i] / count[i];
220 dip[i][XX] = dip[i][XX] / count[i];
221 dip[i][YY] = dip[i][YY] / count[i];
222 dip[i][ZZ] = dip[i][ZZ] / count[i];
226 fprintf(stderr, "No water in slice %d\n", i);
230 *slOrder = sum; /* copy a pointer, I hope */
232 sfree(x0); /* free memory used by coordinate arrays */
235 void h2order_plot(rvec dipole[], real order[], const char *afile,
236 int nslices, real slWidth, const output_env_t oenv)
238 FILE *ord; /* xvgr files with order parameters */
239 int slice; /* loop index */
240 char buf[256]; /* for xvgr title */
241 real factor; /* conversion to Debye from electron*nm */
243 /* factor = 1e-9*1.60217733e-19/3.336e-30 */
244 factor = 1.60217733/3.336e-2;
245 fprintf(stderr, "%d slices\n", nslices);
246 sprintf(buf, "Water orientation with respect to normal");
247 ord = xvgropen(afile, buf,
248 "box (nm)", "mu_x, mu_y, mu_z (D), cosine with normal", oenv);
250 for (slice = 0; slice < nslices; slice++)
252 fprintf(ord, "%8.3f %8.3f %8.3f %8.3f %e\n", slWidth*slice,
253 factor*dipole[slice][XX], factor*dipole[slice][YY],
254 factor*dipole[slice][ZZ], order[slice]);
260 int gmx_h2order(int argc, char *argv[])
262 const char *desc[] = {
263 "[THISMODULE] computes the orientation of water molecules with respect to the normal",
264 "of the box. The program determines the average cosine of the angle",
265 "between the dipole moment of water and an axis of the box. The box is",
266 "divided in slices and the average orientation per slice is printed.",
267 "Each water molecule is assigned to a slice, per time frame, based on the",
268 "position of the oxygen. When [TT]-nm[tt] is used, the angle between the water",
269 "dipole and the axis from the center of mass to the oxygen is calculated",
270 "instead of the angle between the dipole and a box axis."
272 static int axis = 2; /* normal to memb. default z */
273 static const char *axtitle = "Z";
274 static int nslices = 0; /* nr of slices defined */
276 { "-d", FALSE, etSTR, {&axtitle},
277 "Take the normal on the membrane in direction X, Y or Z." },
278 { "-sl", FALSE, etINT, {&nslices},
279 "Calculate order parameter as function of boxlength, dividing the box"
280 " in this number of slices."}
282 const char *bugs[] = {
283 "The program assigns whole water molecules to a slice, based on the first "
284 "atom of three in the index file group. It assumes an order O,H,H. "
285 "Name is not important, but the order is. If this demand is not met, "
286 "assigning molecules to slices is different."
290 real *slOrder, /* av. cosine, per slice */
291 slWidth = 0.0; /* width of a slice */
293 char *grpname, /* groupnames */
295 int ngx, /* nr. of atomsin sol group */
296 nmic = 0; /* nr. of atoms in micelle */
297 t_topology *top; /* topology */
299 atom_id *index, /* indices for solvent group */
301 gmx_bool bMicel = FALSE; /* think we're a micel */
302 t_filenm fnm[] = { /* files for g_order */
303 { efTRX, "-f", NULL, ffREAD }, /* trajectory file */
304 { efNDX, NULL, NULL, ffREAD }, /* index file */
305 { efNDX, "-nm", NULL, ffOPTRD }, /* index with micelle atoms */
306 { efTPR, NULL, NULL, ffREAD }, /* topology file */
307 { efXVG, "-o", "order", ffWRITE }, /* xvgr output file */
310 #define NFILE asize(fnm)
312 if (!parse_common_args(&argc, argv,
313 PCA_CAN_VIEW | PCA_CAN_TIME, NFILE,
314 fnm, asize(pa), pa, asize(desc), desc, asize(bugs), bugs, &oenv))
318 bMicel = opt2bSet("-nm", NFILE, fnm);
320 top = read_top(ftp2fn(efTPR, NFILE, fnm), &ePBC); /* read topology file */
322 rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &ngx, &index, &grpname);
326 rd_index(opt2fn("-nm", NFILE, fnm), 1, &nmic, &micelle, &micname);
329 calc_h2order(ftp2fn(efTRX, NFILE, fnm), index, ngx, &slDipole, &slOrder,
330 &slWidth, &nslices, top, ePBC, axis, bMicel, micelle, nmic,
333 h2order_plot(slDipole, slOrder, opt2fn("-o", NFILE, fnm), nslices,
336 do_view(oenv, opt2fn("-o", NFILE, fnm), "-nxy"); /* view xvgr file */