Merge branch 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 "futil.h"
52 #include "statutil.h"
53 #include "index.h"
54 #include "gmx_ana.h"
55
56
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     */
63 /* directions.                                                              */
64 /****************************************************************************/
65
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)
71 {
72     rvec *x0,            /* coordinates with pbc */
73           dipole,        /* dipole moment due to one molecules */
74           normal,
75           com;           /* center of mass of micel, with bMicel */
76     rvec        *dip;    /* sum of dipoles, unnormalized */
77     matrix       box;    /* box (3x3) */
78     t_trxstatus *status;
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 */
83         i, j, teller = 0,
84         slice = 0,       /* current slice number */
85     *count;              /* nr. of atoms in one slice */
86     gmx_rmpbc_t  gpbc = NULL;
87
88     if ((natoms = read_first_x(oenv, &status, fn, &t, &x0, box)) == 0)
89     {
90         gmx_fatal(FARGS, "Could not read coordinates from statusfile\n");
91     }
92
93     if (!*nslices)
94     {
95         *nslices = (int)(box[axis][axis] * 10); /* default value */
96
97
98     }
99     switch (axis)
100     {
101         case 0:
102             normal[0] = 1; normal[1] = 0; normal[2] = 0;
103             break;
104         case 1:
105             normal[0] = 0; normal[1] = 1; normal[2] = 0;
106             break;
107         case 2:
108             normal[0] = 0; normal[1] = 0; normal[2] = 1;
109             break;
110         default:
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;
114     }
115
116     clear_rvec(dipole);
117     snew(count, *nslices);
118     snew(sum, *nslices);
119     snew(dip, *nslices);
120     snew(frame, *nslices);
121
122     *slWidth = box[axis][axis]/(*nslices);
123     fprintf(stderr, "Box divided in %d slices. Initial width of slice: %f\n",
124             *nslices, *slWidth);
125
126     teller = 0;
127
128     gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms);
129     /*********** Start processing trajectory ***********/
130     do
131     {
132         *slWidth = box[axis][axis]/(*nslices);
133         teller++;
134
135         gmx_rmpbc(gpbc, natoms, box, x0);
136
137         if (bMicel)
138         {
139             calc_xcm(x0, nmic, micel, top->atoms.atom, com, FALSE);
140         }
141
142         for (i = 0; i < ngx/3; i++)
143         {
144             /* put all waters in box */
145             for (j = 0; j < DIM; j++)
146             {
147                 if (x0[index[3*i]][j] < 0)
148                 {
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];
152                 }
153                 if (x0[index[3*i]][j] > box[j][j])
154                 {
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];
158                 }
159             }
160
161             for (j = 0; j < DIM; j++)
162             {
163                 dipole[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;
167             }
168
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.
172              */
173
174             if (bMicel)
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           */
178
179                 sum[slice]   += iprod(dipole, normal) / (norm(dipole) * norm(normal));
180                 frame[slice] += iprod(dipole, normal) / (norm(dipole) * norm(normal));
181                 count[slice]++;
182
183             }
184             else
185             {
186                 /* this is for flat interfaces      */
187
188                 /* determine which slice atom is in */
189                 slice = (x0[index[3*i]][axis] / (*slWidth));
190                 if (slice < 0 || slice >= *nslices)
191                 {
192                     fprintf(stderr, "Coordinate: %f ", x0[index[3*i]][axis]);
193                     fprintf(stderr, "HELP PANIC! slice = %d, OUT OF RANGE!\n", slice);
194                 }
195                 else
196                 {
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 */
202                     count[slice]++;
203                 }
204             }
205         }
206
207     }
208     while (read_next_x(oenv, status, &t, x0, box));
209     /*********** done with status file **********/
210
211     fprintf(stderr, "\nRead trajectory. Printing parameters to file\n");
212     gmx_rmpbc_done(gpbc);
213
214     for (i = 0; i < *nslices; i++) /* average over frames */
215     {
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 */
218         {
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];
223         }
224         else
225         {
226             fprintf(stderr, "No water in slice %d\n", i);
227         }
228     }
229
230     *slOrder  = sum; /* copy a pointer, I hope */
231     *slDipole = dip;
232     sfree(x0);       /* free memory used by coordinate arrays */
233 }
234
235 void h2order_plot(rvec dipole[], real order[], const char *afile,
236                   int nslices, real slWidth, const output_env_t oenv)
237 {
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 */
242
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);
249
250     for (slice = 0; slice < nslices; slice++)
251     {
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]);
255     }
256
257     ffclose(ord);
258 }
259
260 int gmx_h2order(int argc, char *argv[])
261 {
262     const char        *desc[] = {
263         "[TT]g_h2order[tt] 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."
271     };
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       */
275     t_pargs            pa[]    = {
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."}
281     };
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."
287     };
288
289     output_env_t       oenv;
290     real              *slOrder,             /* av. cosine, per slice      */
291                        slWidth = 0.0;       /* width of a slice           */
292     rvec              *slDipole;
293     char              *grpname,             /* groupnames                 */
294     *micname;
295     int                ngx,                 /* nr. of atomsin sol group   */
296                        nmic = 0;            /* nr. of atoms in micelle    */
297     t_topology        *top;                 /* topology           */
298     int                ePBC;
299     atom_id           *index,               /* indices for solvent group  */
300     *micelle                  = NULL;
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         { efTPX, NULL, NULL,  ffREAD },     /* topology file              */
307         { efXVG, "-o",  "order", ffWRITE }, /* xvgr output file       */
308     };
309
310 #define NFILE asize(fnm)
311
312     parse_common_args(&argc, argv,
313                       PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE, NFILE,
314                       fnm, asize(pa), pa, asize(desc), desc, asize(bugs), bugs, &oenv);
315     bMicel = opt2bSet("-nm", NFILE, fnm);
316
317     top = read_top(ftp2fn(efTPX, NFILE, fnm), &ePBC); /* read topology file */
318
319     rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &ngx, &index, &grpname);
320
321     if (bMicel)
322     {
323         rd_index(opt2fn("-nm", NFILE, fnm), 1, &nmic, &micelle, &micname);
324     }
325
326     calc_h2order(ftp2fn(efTRX, NFILE, fnm), index, ngx, &slDipole, &slOrder,
327                  &slWidth, &nslices, top, ePBC, axis, bMicel, micelle, nmic,
328                  oenv);
329
330     h2order_plot(slDipole, slOrder, opt2fn("-o", NFILE, fnm), nslices,
331                  slWidth, oenv);
332
333     do_view(oenv, opt2fn("-o", NFILE, fnm), "-nxy"); /* view xvgr file */
334
335     return 0;
336 }