3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
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.
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.
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.
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.
30 * For more info, check our website at http://www.gromacs.org
33 * Green Red Orange Magenta Azure Cyan Skyblue
51 static void calc_com_pbc(int nrefat, t_topology *top, rvec x[], t_pbc *pbc,
52 atom_id index[], rvec xref, gmx_bool bPBC)
54 const real tol = 1e-4;
60 /* First simple calculation */
63 for (m = 0; (m < nrefat); m++)
66 mass = top->atoms.atom[ai].m;
67 for (j = 0; (j < DIM); j++)
69 xref[j] += mass*x[ai][j];
73 svmul(1/mtot, xref, xref);
74 /* Now check if any atom is more than half the box from the COM */
81 for (m = 0; (m < nrefat); m++)
84 mass = top->atoms.atom[ai].m/mtot;
85 pbc_dx(pbc, x[ai], xref, dx);
86 rvec_add(xref, dx, xtest);
87 for (j = 0; (j < DIM); j++)
89 if (fabs(xtest[j]-x[ai][j]) > tol)
91 /* Here we have used the wrong image for contributing to the COM */
92 xref[j] += mass*(xtest[j]-x[ai][j]);
100 printf("COM: %8.3f %8.3f %8.3f iter = %d\n", xref[XX], xref[YY], xref[ZZ], iter);
108 int gmx_sorient(int argc, char *argv[])
120 int i, j, p, sa0, sa1, sa2, n, ntot, nf, m, *hist1, *hist2, *histn, nbin1, nbin2, nrbin;
121 real *histi1, *histi2, invbw, invrbw;
123 int *isize, nrefgrp, nrefat;
126 real inp, outp, two_pi, nav, normfac, rmin2, rmax2, rcut, rcut2, r2, r, mass, mtot;
130 rvec xref, dx, dxh1, dxh2, outer;
131 gmx_rmpbc_t gpbc = NULL;
133 const char *legr[] = {
134 "<cos(\\8q\\4\\s1\\N)>",
135 "<3cos\\S2\\N(\\8q\\4\\s2\\N)-1>"
137 const char *legc[] = {
138 "cos(\\8q\\4\\s1\\N)",
139 "3cos\\S2\\N(\\8q\\4\\s2\\N)-1"
142 const char *desc[] = {
143 "[TT]g_sorient[tt] analyzes solvent orientation around solutes.",
144 "It calculates two angles between the vector from one or more",
145 "reference positions to the first atom of each solvent molecule:[PAR]",
146 "[GRK]theta[grk][SUB]1[sub]: the angle with the vector from the first atom of the solvent",
147 "molecule to the midpoint between atoms 2 and 3.[BR]",
148 "[GRK]theta[grk][SUB]2[sub]: the angle with the normal of the solvent plane, defined by the",
149 "same three atoms, or, when the option [TT]-v23[tt] is set, ",
150 "the angle with the vector between atoms 2 and 3.[PAR]",
151 "The reference can be a set of atoms or",
152 "the center of mass of a set of atoms. The group of solvent atoms should",
153 "consist of 3 atoms per solvent molecule.",
154 "Only solvent molecules between [TT]-rmin[tt] and [TT]-rmax[tt] are",
155 "considered for [TT]-o[tt] and [TT]-no[tt] each frame.[PAR]",
156 "[TT]-o[tt]: distribtion of [MATH][COS][GRK]theta[grk][SUB]1[sub][cos][math] for rmin<=r<=rmax.[PAR]",
157 "[TT]-no[tt]: distribution of [MATH][COS][GRK]theta[grk][SUB]2[sub][cos][math] for rmin<=r<=rmax.[PAR]",
158 "[TT]-ro[tt]: [MATH][CHEVRON][COS][GRK]theta[grk][SUB]1[sub][cos][chevron][math] and [MATH][CHEVRON]3[COS]^2[GRK]theta[grk][SUB]2[sub][cos]-1[chevron][math] as a function of the",
160 "[TT]-co[tt]: the sum over all solvent molecules within distance r",
161 "of [MATH][COS][GRK]theta[grk][SUB]1[sub][cos][math] and [MATH]3[COS]^2([GRK]theta[grk][SUB]2[sub])-1[cos][math] as a function of r.[PAR]",
162 "[TT]-rc[tt]: the distribution of the solvent molecules as a function of r"
166 static gmx_bool bCom = FALSE, bVec23 = FALSE, bPBC = FALSE;
167 static real rmin = 0.0, rmax = 0.5, binwidth = 0.02, rbinw = 0.02;
169 { "-com", FALSE, etBOOL, {&bCom},
170 "Use the center of mass as the reference postion" },
171 { "-v23", FALSE, etBOOL, {&bVec23},
172 "Use the vector between atoms 2 and 3" },
173 { "-rmin", FALSE, etREAL, {&rmin}, "Minimum distance (nm)" },
174 { "-rmax", FALSE, etREAL, {&rmax}, "Maximum distance (nm)" },
175 { "-cbin", FALSE, etREAL, {&binwidth}, "Binwidth for the cosine" },
176 { "-rbin", FALSE, etREAL, {&rbinw}, "Binwidth for r (nm)" },
177 { "-pbc", FALSE, etBOOL, {&bPBC}, "Check PBC for the center of mass calculation. Only necessary when your reference group consists of several molecules." }
181 { efTRX, NULL, NULL, ffREAD },
182 { efTPS, NULL, NULL, ffREAD },
183 { efNDX, NULL, NULL, ffOPTRD },
184 { efXVG, NULL, "sori.xvg", ffWRITE },
185 { efXVG, "-no", "snor.xvg", ffWRITE },
186 { efXVG, "-ro", "sord.xvg", ffWRITE },
187 { efXVG, "-co", "scum.xvg", ffWRITE },
188 { efXVG, "-rc", "scount.xvg", ffWRITE }
190 #define NFILE asize(fnm)
192 parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW | PCA_BE_NICE,
193 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv);
197 bTPS = (opt2bSet("-s", NFILE, fnm) || !opt2bSet("-n", NFILE, fnm) || bCom);
200 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xtop, NULL, box,
204 /* get index groups */
205 printf("Select a group of reference particles and a solvent group:\n");
211 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 2, isize, index, grpname);
215 get_index(NULL, ftp2fn(efNDX, NFILE, fnm), 2, isize, index, grpname);
231 gmx_fatal(FARGS, "The number of solvent atoms (%d) is not a multiple of 3",
235 /* initialize reading trajectory: */
236 natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
240 rcut = 0.99*sqrt(max_cutoff2(guess_ePBC(box), box));
248 nbin1 = 1+(int)(2*invbw + 0.5);
249 nbin2 = 1+(int)(invbw + 0.5);
255 nrbin = 1+(int)(rcut/rbinw);
271 /* make molecules whole again */
272 gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms);
274 /* start analysis of trajectory */
279 /* make molecules whole again */
280 gmx_rmpbc(gpbc, natoms, box, x);
283 set_pbc(&pbc, ePBC, box);
287 for (p = 0; (p < nrefgrp); p++)
291 calc_com_pbc(nrefat, &top, x, &pbc, index[0], xref, bPBC);
295 copy_rvec(x[index[0][p]], xref);
298 for (m = 0; m < isize[1]; m += 3)
303 range_check(sa0, 0, natoms);
304 range_check(sa1, 0, natoms);
305 range_check(sa2, 0, natoms);
306 pbc_dx(&pbc, x[sa0], xref, dx);
313 /* Determine the normal to the plain */
314 rvec_sub(x[sa1], x[sa0], dxh1);
315 rvec_sub(x[sa2], x[sa0], dxh2);
316 rvec_inc(dxh1, dxh2);
319 inp = iprod(dx, dxh1);
320 cprod(dxh1, dxh2, outer);
322 outp = iprod(dx, outer);
326 /* Use the vector between the 2nd and 3rd atom */
327 rvec_sub(x[sa2], x[sa1], dxh2);
329 outp = iprod(dx, dxh2)/r;
332 int ii = (int)(invrbw*r);
333 range_check(ii, 0, nrbin);
335 histi2[ii] += 3*sqr(outp) - 1;
338 if ((r2 >= rmin2) && (r2 < rmax2))
340 int ii1 = (int)(invbw*(inp + 1));
341 int ii2 = (int)(invbw*fabs(outp));
343 range_check(ii1, 0, nbin1);
344 range_check(ii2, 0, nbin2);
358 while (read_next_x(oenv, status, &t, x, box));
363 gmx_rmpbc_done(gpbc);
365 /* Add the bin for the exact maximum to the previous bin */
366 hist1[nbin1-1] += hist1[nbin1];
367 hist2[nbin2-1] += hist2[nbin2];
369 nav = (real)ntot/(nrefgrp*nf);
370 normfac = invbw/ntot;
372 fprintf(stderr, "Average nr of molecules between %g and %g nm: %.1f\n",
378 fprintf(stderr, "Average cos(theta1) between %g and %g nm: %6.3f\n",
380 fprintf(stderr, "Average 3cos2(theta2)-1 between %g and %g nm: %6.3f\n",
384 sprintf(str, "Solvent orientation between %g and %g nm", rmin, rmax);
385 fp = xvgropen(opt2fn("-o", NFILE, fnm), str, "cos(\\8q\\4\\s1\\N)", "", oenv);
386 if (output_env_get_print_xvgr_codes(oenv))
388 fprintf(fp, "@ subtitle \"average shell size %.1f molecules\"\n", nav);
390 for (i = 0; i < nbin1; i++)
392 fprintf(fp, "%g %g\n", (i+0.5)*binwidth-1, 2*normfac*hist1[i]);
396 sprintf(str, "Solvent normal orientation between %g and %g nm", rmin, rmax);
397 fp = xvgropen(opt2fn("-no", NFILE, fnm), str, "cos(\\8q\\4\\s2\\N)", "", oenv);
398 if (output_env_get_print_xvgr_codes(oenv))
400 fprintf(fp, "@ subtitle \"average shell size %.1f molecules\"\n", nav);
402 for (i = 0; i < nbin2; i++)
404 fprintf(fp, "%g %g\n", (i+0.5)*binwidth, normfac*hist2[i]);
409 sprintf(str, "Solvent orientation");
410 fp = xvgropen(opt2fn("-ro", NFILE, fnm), str, "r (nm)", "", oenv);
411 if (output_env_get_print_xvgr_codes(oenv))
413 fprintf(fp, "@ subtitle \"as a function of distance\"\n");
415 xvgr_legend(fp, 2, legr, oenv);
416 for (i = 0; i < nrbin; i++)
418 fprintf(fp, "%g %g %g\n", (i+0.5)*rbinw,
419 histn[i] ? histi1[i]/histn[i] : 0,
420 histn[i] ? histi2[i]/histn[i] : 0);
424 sprintf(str, "Cumulative solvent orientation");
425 fp = xvgropen(opt2fn("-co", NFILE, fnm), str, "r (nm)", "", oenv);
426 if (output_env_get_print_xvgr_codes(oenv))
428 fprintf(fp, "@ subtitle \"as a function of distance\"\n");
430 xvgr_legend(fp, 2, legc, oenv);
431 normfac = 1.0/(nrefgrp*nf);
434 fprintf(fp, "%g %g %g\n", 0.0, c1, c2);
435 for (i = 0; i < nrbin; i++)
437 c1 += histi1[i]*normfac;
438 c2 += histi2[i]*normfac;
439 fprintf(fp, "%g %g %g\n", (i+1)*rbinw, c1, c2);
443 sprintf(str, "Solvent distribution");
444 fp = xvgropen(opt2fn("-rc", NFILE, fnm), str, "r (nm)", "molecules/nm", oenv);
445 if (output_env_get_print_xvgr_codes(oenv))
447 fprintf(fp, "@ subtitle \"as a function of distance\"\n");
449 normfac = 1.0/(rbinw*nf);
450 for (i = 0; i < nrbin; i++)
452 fprintf(fp, "%g %g\n", (i+0.5)*rbinw, histn[i]*normfac);
456 do_view(oenv, opt2fn("-o", NFILE, fnm), NULL);
457 do_view(oenv, opt2fn("-no", NFILE, fnm), NULL);
458 do_view(oenv, opt2fn("-ro", NFILE, fnm), "-nxy");
459 do_view(oenv, opt2fn("-co", NFILE, fnm), "-nxy");