Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_h2order.c
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.
14
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.
19  *
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.
26  *
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.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Green Red Orange Magenta Azure Cyan Skyblue
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <math.h>
40
41 #include "sysstuff.h"
42 #include <string.h>
43 #include "typedefs.h"
44 #include "smalloc.h"
45 #include "macros.h"
46 #include "princ.h"
47 #include "rmpbc.h"
48 #include "vec.h"
49 #include "xvgr.h"
50 #include "pbc.h"
51 #include "gromacs/fileio/futil.h"
52 #include "statutil.h"
53 #include "index.h"
54 #include "gmx_ana.h"
55 #include "gromacs/fileio/trxio.h"
56
57
58 /****************************************************************************/
59 /* This program calculates the ordering of water molecules across a box, as */
60 /* function of the z-coordinate. This implies averaging over slices and over*/
61 /* time. Output is the average cos of the angle of the dipole with the      */
62 /* normal, per slice.                                                       */
63 /* In addition, it calculates the average dipole moment itself in three     */
64 /* directions.                                                              */
65 /****************************************************************************/
66
67 void calc_h2order(const char *fn, atom_id index[], int ngx, rvec **slDipole,
68                   real **slOrder, real *slWidth, int *nslices,
69                   t_topology *top, int ePBC,
70                   int axis, gmx_bool bMicel, atom_id micel[], int nmic,
71                   const output_env_t oenv)
72 {
73     rvec *x0,            /* coordinates with pbc */
74           dipole,        /* dipole moment due to one molecules */
75           normal,
76           com;           /* center of mass of micel, with bMicel */
77     rvec        *dip;    /* sum of dipoles, unnormalized */
78     matrix       box;    /* box (3x3) */
79     t_trxstatus *status;
80     real         t,      /* time from trajectory */
81     *sum,                /* sum of all cosines of dipoles, per slice */
82     *frame;              /* order over one frame */
83     int natoms,          /* nr. atoms in trj */
84         i, j, teller = 0,
85         slice = 0,       /* current slice number */
86     *count;              /* nr. of atoms in one slice */
87     gmx_rmpbc_t  gpbc = NULL;
88
89     if ((natoms = read_first_x(oenv, &status, fn, &t, &x0, box)) == 0)
90     {
91         gmx_fatal(FARGS, "Could not read coordinates from statusfile\n");
92     }
93
94     if (!*nslices)
95     {
96         *nslices = (int)(box[axis][axis] * 10); /* default value */
97
98
99     }
100     switch (axis)
101     {
102         case 0:
103             normal[0] = 1; normal[1] = 0; normal[2] = 0;
104             break;
105         case 1:
106             normal[0] = 0; normal[1] = 1; normal[2] = 0;
107             break;
108         case 2:
109             normal[0] = 0; normal[1] = 0; normal[2] = 1;
110             break;
111         default:
112             gmx_fatal(FARGS, "No valid value for -axis-. Exiting.\n");
113             /* make compiler happy */
114             normal[0] = 1; normal[1] = 0; normal[2] = 0;
115     }
116
117     clear_rvec(dipole);
118     snew(count, *nslices);
119     snew(sum, *nslices);
120     snew(dip, *nslices);
121     snew(frame, *nslices);
122
123     *slWidth = box[axis][axis]/(*nslices);
124     fprintf(stderr, "Box divided in %d slices. Initial width of slice: %f\n",
125             *nslices, *slWidth);
126
127     teller = 0;
128
129     gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms);
130     /*********** Start processing trajectory ***********/
131     do
132     {
133         *slWidth = box[axis][axis]/(*nslices);
134         teller++;
135
136         gmx_rmpbc(gpbc, natoms, box, x0);
137
138         if (bMicel)
139         {
140             calc_xcm(x0, nmic, micel, top->atoms.atom, com, FALSE);
141         }
142
143         for (i = 0; i < ngx/3; i++)
144         {
145             /* put all waters in box */
146             for (j = 0; j < DIM; j++)
147             {
148                 if (x0[index[3*i]][j] < 0)
149                 {
150                     x0[index[3*i]][j]   += box[j][j];
151                     x0[index[3*i+1]][j] += box[j][j];
152                     x0[index[3*i+2]][j] += box[j][j];
153                 }
154                 if (x0[index[3*i]][j] > box[j][j])
155                 {
156                     x0[index[3*i]][j]   -= box[j][j];
157                     x0[index[3*i+1]][j] -= box[j][j];
158                     x0[index[3*i+2]][j] -= box[j][j];
159                 }
160             }
161
162             for (j = 0; j < DIM; j++)
163             {
164                 dipole[j] =
165                     x0[index[3*i]][j] * top->atoms.atom[index[3*i]].q +
166                     x0[index[3*i+1]][j] * top->atoms.atom[index[3*i+1]].q +
167                     x0[index[3*i+2]][j] * top->atoms.atom[index[3*i+2]].q;
168             }
169
170             /* now we have a dipole vector. Might as well safe it. Then the
171                rest depends on whether we're dealing with
172                a flat or a spherical interface.
173              */
174
175             if (bMicel)
176             {                                          /* this is for spherical interfaces */
177                 rvec_sub(com, x0[index[3*i]], normal); /* vector from Oxygen to COM */
178                 slice = norm(normal)/(*slWidth);       /* spherical slice           */
179
180                 sum[slice]   += iprod(dipole, normal) / (norm(dipole) * norm(normal));
181                 frame[slice] += iprod(dipole, normal) / (norm(dipole) * norm(normal));
182                 count[slice]++;
183
184             }
185             else
186             {
187                 /* this is for flat interfaces      */
188
189                 /* determine which slice atom is in */
190                 slice = (x0[index[3*i]][axis] / (*slWidth));
191                 if (slice < 0 || slice >= *nslices)
192                 {
193                     fprintf(stderr, "Coordinate: %f ", x0[index[3*i]][axis]);
194                     fprintf(stderr, "HELP PANIC! slice = %d, OUT OF RANGE!\n", slice);
195                 }
196                 else
197                 {
198                     rvec_add(dipole, dip[slice], dip[slice]);
199                     /* Add dipole to total. mag[slice] is total dipole in axis direction */
200                     sum[slice]   += iprod(dipole, normal)/norm(dipole);
201                     frame[slice] += iprod(dipole, normal)/norm(dipole);
202                     /* increase count for that slice */
203                     count[slice]++;
204                 }
205             }
206         }
207
208     }
209     while (read_next_x(oenv, status, &t, x0, box));
210     /*********** done with status file **********/
211
212     fprintf(stderr, "\nRead trajectory. Printing parameters to file\n");
213     gmx_rmpbc_done(gpbc);
214
215     for (i = 0; i < *nslices; i++) /* average over frames */
216     {
217         fprintf(stderr, "%d waters in slice %d\n", count[i], i);
218         if (count[i] > 0) /* divide by number of molecules in each slice */
219         {
220             sum[i]     = sum[i] / count[i];
221             dip[i][XX] = dip[i][XX] / count[i];
222             dip[i][YY] = dip[i][YY] / count[i];
223             dip[i][ZZ] = dip[i][ZZ] / count[i];
224         }
225         else
226         {
227             fprintf(stderr, "No water in slice %d\n", i);
228         }
229     }
230
231     *slOrder  = sum; /* copy a pointer, I hope */
232     *slDipole = dip;
233     sfree(x0);       /* free memory used by coordinate arrays */
234 }
235
236 void h2order_plot(rvec dipole[], real order[], const char *afile,
237                   int nslices, real slWidth, const output_env_t oenv)
238 {
239     FILE       *ord;              /* xvgr files with order parameters  */
240     int         slice;            /* loop index     */
241     char        buf[256];         /* for xvgr title */
242     real        factor;           /* conversion to Debye from electron*nm */
243
244     /*  factor = 1e-9*1.60217733e-19/3.336e-30 */
245     factor = 1.60217733/3.336e-2;
246     fprintf(stderr, "%d slices\n", nslices);
247     sprintf(buf, "Water orientation with respect to normal");
248     ord = xvgropen(afile, buf,
249                    "box (nm)", "mu_x, mu_y, mu_z (D), cosine with normal", oenv);
250
251     for (slice = 0; slice < nslices; slice++)
252     {
253         fprintf(ord, "%8.3f %8.3f %8.3f %8.3f %e\n", slWidth*slice,
254                 factor*dipole[slice][XX], factor*dipole[slice][YY],
255                 factor*dipole[slice][ZZ], order[slice]);
256     }
257
258     ffclose(ord);
259 }
260
261 int gmx_h2order(int argc, char *argv[])
262 {
263     const char        *desc[] = {
264         "[TT]g_h2order[tt] computes the orientation of water molecules with respect to the normal",
265         "of the box. The program determines the average cosine of the angle",
266         "between the dipole moment of water and an axis of the box. The box is",
267         "divided in slices and the average orientation per slice is printed.",
268         "Each water molecule is assigned to a slice, per time frame, based on the",
269         "position of the oxygen. When [TT]-nm[tt] is used, the angle between the water",
270         "dipole and the axis from the center of mass to the oxygen is calculated",
271         "instead of the angle between the dipole and a box axis."
272     };
273     static int         axis    = 2;           /* normal to memb. default z  */
274     static const char *axtitle = "Z";
275     static int         nslices = 0;           /* nr of slices defined       */
276     t_pargs            pa[]    = {
277         { "-d",   FALSE, etSTR, {&axtitle},
278           "Take the normal on the membrane in direction X, Y or Z." },
279         { "-sl",  FALSE, etINT, {&nslices},
280           "Calculate order parameter as function of boxlength, dividing the box"
281           " in this number of slices."}
282     };
283     const char        *bugs[] = {
284         "The program assigns whole water molecules to a slice, based on the first "
285         "atom of three in the index file group. It assumes an order O,H,H. "
286         "Name is not important, but the order is. If this demand is not met, "
287         "assigning molecules to slices is different."
288     };
289
290     output_env_t       oenv;
291     real              *slOrder,             /* av. cosine, per slice      */
292                        slWidth = 0.0;       /* width of a slice           */
293     rvec              *slDipole;
294     char              *grpname,             /* groupnames                 */
295     *micname;
296     int                ngx,                 /* nr. of atomsin sol group   */
297                        nmic = 0;            /* nr. of atoms in micelle    */
298     t_topology        *top;                 /* topology           */
299     int                ePBC;
300     atom_id           *index,               /* indices for solvent group  */
301     *micelle                  = NULL;
302     gmx_bool           bMicel =  FALSE;     /* think we're a micel        */
303     t_filenm           fnm[]  = {           /* files for g_order      */
304         { efTRX, "-f", NULL,  ffREAD },     /* trajectory file            */
305         { efNDX, NULL, NULL,  ffREAD },     /* index file         */
306         { efNDX, "-nm", NULL, ffOPTRD },    /* index with micelle atoms   */
307         { efTPX, NULL, NULL,  ffREAD },     /* topology file              */
308         { efXVG, "-o",  "order", ffWRITE }, /* xvgr output file       */
309     };
310
311 #define NFILE asize(fnm)
312
313     if (!parse_common_args(&argc, argv,
314                            PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE, NFILE,
315                            fnm, asize(pa), pa, asize(desc), desc, asize(bugs), bugs, &oenv))
316     {
317         return 0;
318     }
319     bMicel = opt2bSet("-nm", NFILE, fnm);
320
321     top = read_top(ftp2fn(efTPX, NFILE, fnm), &ePBC); /* read topology file */
322
323     rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &ngx, &index, &grpname);
324
325     if (bMicel)
326     {
327         rd_index(opt2fn("-nm", NFILE, fnm), 1, &nmic, &micelle, &micname);
328     }
329
330     calc_h2order(ftp2fn(efTRX, NFILE, fnm), index, ngx, &slDipole, &slOrder,
331                  &slWidth, &nslices, top, ePBC, axis, bMicel, micelle, nmic,
332                  oenv);
333
334     h2order_plot(slDipole, slOrder, opt2fn("-o", NFILE, fnm), nslices,
335                  slWidth, oenv);
336
337     do_view(oenv, opt2fn("-o", NFILE, fnm), "-nxy"); /* view xvgr file */
338
339     return 0;
340 }