Remove gmx custom fixed int (e.g. gmx_int64_t) types
[gromacs.git] / src / gromacs / gmxana / gmx_h2order.cpp
blobc5edc450a65d9e972f0ad95f4c7f3d0935468a99
1 /*
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,2015,2016,2017,2018, 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.
37 #include "gmxpre.h"
39 #include <cmath>
40 #include <cstring>
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"
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 /****************************************************************************/
66 static void calc_h2order(const char *fn, const int index[], int ngx, rvec **slDipole,
67 real **slOrder, real *slWidth, int *nslices,
68 const t_topology *top, int ePBC,
69 int axis, gmx_bool bMicel, int micel[], int nmic,
70 const gmx_output_env_t *oenv)
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 = nullptr;
88 if ((natoms = read_first_x(oenv, &status, fn, &t, &x0, box)) == 0)
90 gmx_fatal(FARGS, "Could not read coordinates from statusfile\n");
93 if (!*nslices)
95 *nslices = static_cast<int>(box[axis][axis] * 10); /* default value */
99 switch (axis)
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");
114 clear_rvec(dipole);
115 snew(count, *nslices);
116 snew(sum, *nslices);
117 snew(dip, *nslices);
118 snew(frame, *nslices);
120 *slWidth = box[axis][axis]/(*nslices);
121 fprintf(stderr, "Box divided in %d slices. Initial width of slice: %f\n",
122 *nslices, *slWidth);
124 teller = 0;
126 gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms);
127 /*********** Start processing trajectory ***********/
130 *slWidth = box[axis][axis]/(*nslices);
131 teller++;
133 gmx_rmpbc(gpbc, natoms, box, x0);
135 if (bMicel)
137 calc_xcm(x0, nmic, micel, top->atoms.atom, com, FALSE);
140 for (i = 0; i < ngx/3; i++)
142 /* put all waters in box */
143 for (j = 0; j < DIM; j++)
145 if (x0[index[3*i]][j] < 0)
147 x0[index[3*i]][j] += box[j][j];
148 x0[index[3*i+1]][j] += box[j][j];
149 x0[index[3*i+2]][j] += box[j][j];
151 if (x0[index[3*i]][j] > box[j][j])
153 x0[index[3*i]][j] -= box[j][j];
154 x0[index[3*i+1]][j] -= box[j][j];
155 x0[index[3*i+2]][j] -= box[j][j];
159 for (j = 0; j < DIM; j++)
161 dipole[j] =
162 x0[index[3*i]][j] * top->atoms.atom[index[3*i]].q +
163 x0[index[3*i+1]][j] * top->atoms.atom[index[3*i+1]].q +
164 x0[index[3*i+2]][j] * top->atoms.atom[index[3*i+2]].q;
167 /* now we have a dipole vector. Might as well safe it. Then the
168 rest depends on whether we're dealing with
169 a flat or a spherical interface.
172 if (bMicel)
173 { /* this is for spherical interfaces */
174 rvec_sub(com, x0[index[3*i]], normal); /* vector from Oxygen to COM */
175 slice = static_cast<int>(norm(normal)/(*slWidth)); /* spherical slice */
177 sum[slice] += iprod(dipole, normal) / (norm(dipole) * norm(normal));
178 frame[slice] += iprod(dipole, normal) / (norm(dipole) * norm(normal));
179 count[slice]++;
182 else
184 /* this is for flat interfaces */
186 /* determine which slice atom is in */
187 slice = static_cast<int>(x0[index[3*i]][axis] / (*slWidth));
188 if (slice < 0 || slice >= *nslices)
190 fprintf(stderr, "Coordinate: %f ", x0[index[3*i]][axis]);
191 fprintf(stderr, "HELP PANIC! slice = %d, OUT OF RANGE!\n", slice);
193 else
195 rvec_add(dipole, dip[slice], dip[slice]);
196 /* Add dipole to total. mag[slice] is total dipole in axis direction */
197 sum[slice] += iprod(dipole, normal)/norm(dipole);
198 frame[slice] += iprod(dipole, normal)/norm(dipole);
199 /* increase count for that slice */
200 count[slice]++;
206 while (read_next_x(oenv, status, &t, x0, box));
207 /*********** done with status file **********/
209 fprintf(stderr, "\nRead trajectory. Printing parameters to file\n");
210 gmx_rmpbc_done(gpbc);
212 for (i = 0; i < *nslices; i++) /* average over frames */
214 fprintf(stderr, "%d waters in slice %d\n", count[i], i);
215 if (count[i] > 0) /* divide by number of molecules in each slice */
217 sum[i] = sum[i] / count[i];
218 dip[i][XX] = dip[i][XX] / count[i];
219 dip[i][YY] = dip[i][YY] / count[i];
220 dip[i][ZZ] = dip[i][ZZ] / count[i];
222 else
224 fprintf(stderr, "No water in slice %d\n", i);
228 *slOrder = sum; /* copy a pointer, I hope */
229 *slDipole = dip;
230 sfree(x0); /* free memory used by coordinate arrays */
233 static void h2order_plot(rvec dipole[], real order[], const char *afile,
234 int nslices, real slWidth, const gmx_output_env_t *oenv)
236 FILE *ord; /* xvgr files with order parameters */
237 int slice; /* loop index */
238 char buf[256]; /* for xvgr title */
239 real factor; /* conversion to Debye from electron*nm */
241 /* factor = 1e-9*1.60217733e-19/3.336e-30 */
242 factor = 1.60217733/3.336e-2;
243 fprintf(stderr, "%d slices\n", nslices);
244 sprintf(buf, "Water orientation with respect to normal");
245 ord = xvgropen(afile, buf,
246 "box (nm)", "mu_x, mu_y, mu_z (D), cosine with normal", oenv);
248 for (slice = 0; slice < nslices; slice++)
250 fprintf(ord, "%8.3f %8.3f %8.3f %8.3f %e\n", slWidth*slice,
251 factor*dipole[slice][XX], factor*dipole[slice][YY],
252 factor*dipole[slice][ZZ], order[slice]);
255 xvgrclose(ord);
258 int gmx_h2order(int argc, char *argv[])
260 const char *desc[] = {
261 "[THISMODULE] computes the orientation of water molecules with respect to the normal",
262 "of the box. The program determines the average cosine of the angle",
263 "between the dipole moment of water and an axis of the box. The box is",
264 "divided in slices and the average orientation per slice is printed.",
265 "Each water molecule is assigned to a slice, per time frame, based on the",
266 "position of the oxygen. When [TT]-nm[tt] is used, the angle between the water",
267 "dipole and the axis from the center of mass to the oxygen is calculated",
268 "instead of the angle between the dipole and a box axis."
270 static int axis = 2; /* normal to memb. default z */
271 static const char *axtitle = "Z";
272 static int nslices = 0; /* nr of slices defined */
273 t_pargs pa[] = {
274 { "-d", FALSE, etSTR, {&axtitle},
275 "Take the normal on the membrane in direction X, Y or Z." },
276 { "-sl", FALSE, etINT, {&nslices},
277 "Calculate order parameter as function of boxlength, dividing the box"
278 " in this number of slices."}
280 const char *bugs[] = {
281 "The program assigns whole water molecules to a slice, based on the first "
282 "atom of three in the index file group. It assumes an order O,H,H. "
283 "Name is not important, but the order is. If this demand is not met, "
284 "assigning molecules to slices is different."
287 gmx_output_env_t *oenv;
288 real *slOrder, /* av. cosine, per slice */
289 slWidth = 0.0; /* width of a slice */
290 rvec *slDipole;
291 char *grpname, /* groupnames */
292 *micname;
293 int ngx, /* nr. of atomsin sol group */
294 nmic = 0; /* nr. of atoms in micelle */
295 t_topology *top; /* topology */
296 int ePBC;
297 int *index, /* indices for solvent group */
298 *micelle = nullptr;
299 gmx_bool bMicel = FALSE; /* think we're a micel */
300 t_filenm fnm[] = { /* files for g_order */
301 { efTRX, "-f", nullptr, ffREAD }, /* trajectory file */
302 { efNDX, nullptr, nullptr, ffREAD }, /* index file */
303 { efNDX, "-nm", nullptr, ffOPTRD }, /* index with micelle atoms */
304 { efTPR, nullptr, nullptr, ffREAD }, /* topology file */
305 { efXVG, "-o", "order", ffWRITE }, /* xvgr output file */
308 #define NFILE asize(fnm)
310 if (!parse_common_args(&argc, argv,
311 PCA_CAN_VIEW | PCA_CAN_TIME, NFILE,
312 fnm, asize(pa), pa, asize(desc), desc, asize(bugs), bugs, &oenv))
314 return 0;
316 bMicel = opt2bSet("-nm", NFILE, fnm);
318 top = read_top(ftp2fn(efTPR, NFILE, fnm), &ePBC); /* read topology file */
320 rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &ngx, &index, &grpname);
322 if (bMicel)
324 rd_index(opt2fn("-nm", NFILE, fnm), 1, &nmic, &micelle, &micname);
327 calc_h2order(ftp2fn(efTRX, NFILE, fnm), index, ngx, &slDipole, &slOrder,
328 &slWidth, &nslices, top, ePBC, axis, bMicel, micelle, nmic,
329 oenv);
331 h2order_plot(slDipole, slOrder, opt2fn("-o", NFILE, fnm), nslices,
332 slWidth, oenv);
334 do_view(oenv, opt2fn("-o", NFILE, fnm), "-nxy"); /* view xvgr file */
336 return 0;