Fix gmx h2order -d option
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_h2order.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2013,2014,2015,2016,2017,2018,2019,2020, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include <cmath>
40 #include <cstring>
41
42 #include "gromacs/commandline/pargs.h"
43 #include "gromacs/commandline/viewit.h"
44 #include "gromacs/fileio/trxio.h"
45 #include "gromacs/fileio/xvgr.h"
46 #include "gromacs/gmxana/gmx_ana.h"
47 #include "gromacs/gmxana/princ.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/pbcutil/rmpbc.h"
50 #include "gromacs/topology/index.h"
51 #include "gromacs/topology/topology.h"
52 #include "gromacs/utility/arraysize.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/futil.h"
55 #include "gromacs/utility/smalloc.h"
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 static void calc_h2order(const char*             fn,
67                          const int               index[],
68                          int                     ngx,
69                          rvec**                  slDipole,
70                          real**                  slOrder,
71                          real*                   slWidth,
72                          int*                    nslices,
73                          const t_topology*       top,
74                          int                     ePBC,
75                          int                     axis,
76                          gmx_bool                bMicel,
77                          int                     micel[],
78                          int                     nmic,
79                          const gmx_output_env_t* oenv)
80 {
81     rvec *x0,            /* coordinates with pbc */
82             dipole,      /* dipole moment due to one molecules */
83             normal, com; /* center of mass of micel, with bMicel */
84     rvec*        dip;    /* sum of dipoles, unnormalized */
85     matrix       box;    /* box (3x3) */
86     t_trxstatus* status;
87     real         t,                      /* time from trajectory */
88             *sum,                        /* sum of all cosines of dipoles, per slice */
89             *frame;                      /* order over one frame */
90     int natoms,                          /* nr. atoms in trj */
91             i, j, teller = 0, slice = 0, /* current slice number */
92             *count;                      /* nr. of atoms in one slice */
93     gmx_rmpbc_t gpbc = nullptr;
94
95     if ((natoms = read_first_x(oenv, &status, fn, &t, &x0, box)) == 0)
96     {
97         gmx_fatal(FARGS, "Could not read coordinates from statusfile\n");
98     }
99
100     if (!*nslices)
101     {
102         *nslices = static_cast<int>(box[axis][axis] * 10); /* default value */
103     }
104     switch (axis)
105     {
106         case 0:
107             normal[0] = 1;
108             normal[1] = 0;
109             normal[2] = 0;
110             break;
111         case 1:
112             normal[0] = 0;
113             normal[1] = 1;
114             normal[2] = 0;
115             break;
116         case 2:
117             normal[0] = 0;
118             normal[1] = 0;
119             normal[2] = 1;
120             break;
121         default: gmx_fatal(FARGS, "No valid value for -axis-. Exiting.\n");
122     }
123
124     clear_rvec(dipole);
125     snew(count, *nslices);
126     snew(sum, *nslices);
127     snew(dip, *nslices);
128     snew(frame, *nslices);
129
130     *slWidth = box[axis][axis] / (*nslices);
131     fprintf(stderr, "Box divided in %d slices. Initial width of slice: %f\n", *nslices, *slWidth);
132
133     teller = 0;
134
135     gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms);
136     /*********** Start processing trajectory ***********/
137     do
138     {
139         *slWidth = box[axis][axis] / (*nslices);
140         teller++;
141
142         gmx_rmpbc(gpbc, natoms, box, x0);
143
144         if (bMicel)
145         {
146             calc_xcm(x0, nmic, micel, top->atoms.atom, com, FALSE);
147         }
148
149         for (i = 0; i < ngx / 3; i++)
150         {
151             /* put all waters in box */
152             for (j = 0; j < DIM; j++)
153             {
154                 if (x0[index[3 * i]][j] < 0)
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                 if (x0[index[3 * i]][j] > box[j][j])
161                 {
162                     x0[index[3 * i]][j] -= box[j][j];
163                     x0[index[3 * i + 1]][j] -= box[j][j];
164                     x0[index[3 * i + 2]][j] -= box[j][j];
165                 }
166             }
167
168             for (j = 0; j < DIM; j++)
169             {
170                 dipole[j] = x0[index[3 * i]][j] * top->atoms.atom[index[3 * i]].q
171                             + x0[index[3 * i + 1]][j] * top->atoms.atom[index[3 * i + 1]].q
172                             + x0[index[3 * i + 2]][j] * top->atoms.atom[index[3 * i + 2]].q;
173             }
174
175             /* now we have a dipole vector. Might as well safe it. Then the
176                rest depends on whether we're dealing with
177                a flat or a spherical interface.
178              */
179
180             if (bMicel)
181             {                                            /* this is for spherical interfaces */
182                 rvec_sub(com, x0[index[3 * i]], normal); /* vector from Oxygen to COM */
183                 slice = static_cast<int>(norm(normal) / (*slWidth)); /* spherical slice           */
184
185                 sum[slice] += iprod(dipole, normal) / (norm(dipole) * norm(normal));
186                 frame[slice] += iprod(dipole, normal) / (norm(dipole) * norm(normal));
187                 count[slice]++;
188             }
189             else
190             {
191                 /* this is for flat interfaces      */
192
193                 /* determine which slice atom is in */
194                 slice = static_cast<int>(x0[index[3 * i]][axis] / (*slWidth));
195                 if (slice < 0 || slice >= *nslices)
196                 {
197                     fprintf(stderr, "Coordinate: %f ", x0[index[3 * i]][axis]);
198                     fprintf(stderr, "HELP PANIC! slice = %d, OUT OF RANGE!\n", slice);
199                 }
200                 else
201                 {
202                     rvec_add(dipole, dip[slice], dip[slice]);
203                     /* Add dipole to total. mag[slice] is total dipole in axis direction */
204                     sum[slice] += iprod(dipole, normal) / norm(dipole);
205                     frame[slice] += iprod(dipole, normal) / norm(dipole);
206                     /* increase count for that slice */
207                     count[slice]++;
208                 }
209             }
210         }
211
212     } while (read_next_x(oenv, status, &t, x0, box));
213     /*********** done with status file **********/
214
215     fprintf(stderr, "\nRead trajectory. Printing parameters to file\n");
216     gmx_rmpbc_done(gpbc);
217
218     for (i = 0; i < *nslices; i++) /* average over frames */
219     {
220         fprintf(stderr, "%d waters in slice %d\n", count[i], i);
221         if (count[i] > 0) /* divide by number of molecules in each slice */
222         {
223             sum[i]     = sum[i] / count[i];
224             dip[i][XX] = dip[i][XX] / count[i];
225             dip[i][YY] = dip[i][YY] / count[i];
226             dip[i][ZZ] = dip[i][ZZ] / count[i];
227         }
228         else
229         {
230             fprintf(stderr, "No water in slice %d\n", i);
231         }
232     }
233
234     *slOrder  = sum; /* copy a pointer, I hope */
235     *slDipole = dip;
236     sfree(x0); /* free memory used by coordinate arrays */
237 }
238
239 static void h2order_plot(rvec dipole[], real order[], const char* afile, int nslices, real slWidth, const gmx_output_env_t* oenv)
240 {
241     FILE* ord;      /* xvgr files with order parameters  */
242     int   slice;    /* loop index     */
243     char  buf[256]; /* for xvgr title */
244     real  factor;   /* conversion to Debye from electron*nm */
245
246     /*  factor = 1e-9*1.60217733e-19/3.336e-30 */
247     factor = 1.60217733 / 3.336e-2;
248     fprintf(stderr, "%d slices\n", nslices);
249     sprintf(buf, "Water orientation with respect to normal");
250     ord = xvgropen(afile, buf, "box (nm)", "mu_x, mu_y, mu_z (D), cosine with normal", oenv);
251
252     for (slice = 0; slice < nslices; slice++)
253     {
254         fprintf(ord, "%8.3f %8.3f %8.3f %8.3f %e\n", slWidth * slice, factor * dipole[slice][XX],
255                 factor * dipole[slice][YY], factor * dipole[slice][ZZ], order[slice]);
256     }
257
258     xvgrclose(ord);
259 }
260
261 enum
262 {
263     axisSEL,
264     axisZ,
265     axisY,
266     axisX,
267     axisNR
268 };
269
270 int gmx_h2order(int argc, char* argv[])
271 {
272     const char* desc[] = {
273         "[THISMODULE] computes the orientation of water molecules with respect to the normal",
274         "of the box. The program determines the average cosine of the angle",
275         "between the dipole moment of water and an axis of the box. The box is",
276         "divided in slices and the average orientation per slice is printed.",
277         "Each water molecule is assigned to a slice, per time frame, based on the",
278         "position of the oxygen. When [TT]-nm[tt] is used, the angle between the water",
279         "dipole and the axis from the center of mass to the oxygen is calculated",
280         "instead of the angle between the dipole and a box axis."
281     };
282     static const char* axisOption[axisNR + 1] = { nullptr, "Z", "Y", "X", nullptr };
283     static int         nslices                = 0; /* nr of slices defined       */
284     // The struct that will hold the parsed user input
285     t_pargs     pa[]   = { { "-d",
286                        FALSE,
287                        etENUM,
288                        { axisOption },
289                        "Take the normal on the membrane in direction X, Y or Z." },
290                      { "-sl",
291                        FALSE,
292                        etINT,
293                        { &nslices },
294                        "Calculate order parameter as function of boxlength, dividing the box"
295                        " in this number of slices." } };
296     const char* bugs[] = {
297         "The program assigns whole water molecules to a slice, based on the first "
298         "atom of three in the index file group. It assumes an order O,H,H. "
299         "Name is not important, but the order is. If this demand is not met, "
300         "assigning molecules to slices is different."
301     };
302
303     gmx_output_env_t* oenv;
304     real *            slOrder, /* av. cosine, per slice      */
305             slWidth = 0.0;     /* width of a slice           */
306     rvec* slDipole;
307     char *grpname, /* groupnames                 */
308             *micname;
309     int ngx,          /* nr. of atomsin sol group   */
310             nmic = 0; /* nr. of atoms in micelle    */
311     t_topology* top;  /* topology           */
312     int         ePBC;
313     int *       index, /* indices for solvent group  */
314             *micelle = nullptr;
315     gmx_bool bMicel  = FALSE; /* think we're a micel        */
316     t_filenm fnm[]   = {
317         /* files for g_order      */
318         { efTRX, "-f", nullptr, ffREAD },    /* trajectory file            */
319         { efNDX, nullptr, nullptr, ffREAD }, /* index file         */
320         { efNDX, "-nm", nullptr, ffOPTRD },  /* index with micelle atoms   */
321         { efTPR, nullptr, nullptr, ffREAD }, /* topology file              */
322         { efXVG, "-o", "order", ffWRITE },   /* xvgr output file       */
323     };
324
325 #define NFILE asize(fnm)
326
327     // Parse the user input in argv into pa
328     if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME, NFILE, fnm, asize(pa), pa,
329                            asize(desc), desc, asize(bugs), bugs, &oenv))
330     {
331         return 0;
332     }
333
334     // Process the axis option chosen by the user to set the
335     // axis used for the computation. The useful choice is an
336     // axis normal to the membrane. Default is z-axis.
337     int axis = ZZ;
338     switch (nenum(axisOption))
339     {
340         case axisX: axis = XX; break;
341         case axisY: axis = YY; break;
342         case axisZ: axis = ZZ; break;
343         default: axis = ZZ;
344     }
345
346     bMicel = opt2bSet("-nm", NFILE, fnm);
347
348     top = read_top(ftp2fn(efTPR, NFILE, fnm), &ePBC); /* read topology file */
349
350     rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &ngx, &index, &grpname);
351
352     if (bMicel)
353     {
354         rd_index(opt2fn("-nm", NFILE, fnm), 1, &nmic, &micelle, &micname);
355     }
356
357     calc_h2order(ftp2fn(efTRX, NFILE, fnm), index, ngx, &slDipole, &slOrder, &slWidth, &nslices,
358                  top, ePBC, axis, bMicel, micelle, nmic, oenv);
359
360     h2order_plot(slDipole, slOrder, opt2fn("-o", NFILE, fnm), nslices, slWidth, oenv);
361
362     do_view(oenv, opt2fn("-o", NFILE, fnm), "-nxy"); /* view xvgr file */
363
364     return 0;
365 }