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
62 real calc_gyro(rvec x[], int gnx, atom_id index[], t_atom atom[], real tm,
63 rvec gvec, rvec d, gmx_bool bQ, gmx_bool bRot, gmx_bool bMOI, matrix trans)
66 real gyro, dx2, m0, Itot;
71 principal_comp(gnx, index, atom, x, trans, d);
77 for (m = 0; (m < DIM); m++)
82 pr_rvecs(stderr, 0, "trans", trans, DIM);
84 /* rotate_atoms(gnx,index,x,trans); */
87 for (i = 0; (i < gnx); i++)
92 m0 = fabs(atom[ii].q);
98 for (m = 0; (m < DIM); m++)
100 dx2 = x[ii][m]*x[ii][m];
104 gyro = comp[XX]+comp[YY]+comp[ZZ];
106 for (m = 0; (m < DIM); m++)
108 gvec[m] = sqrt((gyro-comp[m])/tm);
111 return sqrt(gyro/tm);
114 void calc_gyro_z(rvec x[], matrix box,
115 int gnx, atom_id index[], t_atom atom[],
116 int nz, real time, FILE *out)
118 static dvec *inertia = NULL;
119 static double *tm = NULL;
121 real zf, w, sdet, e1, e2;
129 for (i = 0; i < nz; i++)
131 clear_dvec(inertia[i]);
135 for (i = 0; (i < gnx); i++)
138 zf = nz*x[ii][ZZ]/box[ZZ][ZZ];
147 for (j = 0; j < 2; j++)
154 w = atom[ii].m*(1 + cos(M_PI*(zf - zi)));
155 inertia[zi][0] += w*sqr(x[ii][YY]);
156 inertia[zi][1] += w*sqr(x[ii][XX]);
157 inertia[zi][2] -= w*x[ii][XX]*x[ii][YY];
161 fprintf(out, "%10g", time);
162 for (j = 0; j < nz; j++)
164 for (i = 0; i < 3; i++)
166 inertia[j][i] /= tm[j];
168 sdet = sqrt(sqr(inertia[j][0] - inertia[j][1]) + 4*sqr(inertia[j][2]));
169 e1 = sqrt(0.5*(inertia[j][0] + inertia[j][1] + sdet));
170 e2 = sqrt(0.5*(inertia[j][0] + inertia[j][1] - sdet));
171 fprintf(out, " %5.3f %5.3f", e1, e2);
176 int gmx_gyrate(int argc, char *argv[])
178 const char *desc[] = {
179 "[TT]g_gyrate[tt] computes the radius of gyration of a molecule",
180 "and the radii of gyration about the [IT]x[it]-, [IT]y[it]- and [IT]z[it]-axes,",
181 "as a function of time. The atoms are explicitly mass weighted.[PAR]",
182 "With the [TT]-nmol[tt] option the radius of gyration will be calculated",
183 "for multiple molecules by splitting the analysis group in equally",
185 "With the option [TT]-nz[tt] 2D radii of gyration in the [IT]x-y[it] plane",
186 "of slices along the [IT]z[it]-axis are calculated."
188 static int nmol = 1, nz = 0;
189 static gmx_bool bQ = FALSE, bRot = FALSE, bMOI = FALSE;
191 { "-nmol", FALSE, etINT, {&nmol},
192 "The number of molecules to analyze" },
193 { "-q", FALSE, etBOOL, {&bQ},
194 "Use absolute value of the charge of an atom as weighting factor instead of mass" },
195 { "-p", FALSE, etBOOL, {&bRot},
196 "Calculate the radii of gyration about the principal axes." },
197 { "-moi", FALSE, etBOOL, {&bMOI},
198 "Calculate the moments of inertia (defined by the principal axes)." },
199 { "-nz", FALSE, etINT, {&nz},
200 "Calculate the 2D radii of gyration of this number of slices along the z-axis" },
207 rvec xcm, gvec, gvec1;
210 real **moi_trans = NULL;
211 int max_moi = 0, delta_moi = 100;
212 rvec d, d1; /* eigenvalues of inertia tensor */
213 real t, t0, tm, gyro;
215 char *grpname, title[256];
216 int i, j, m, gnx, nam, mol;
219 gmx_rmpbc_t gpbc = NULL;
220 const char *leg[] = { "Rg", "RgX", "RgY", "RgZ" };
221 const char *legI[] = { "Itot", "I1", "I2", "I3" };
222 #define NLEG asize(leg)
224 { efTRX, "-f", NULL, ffREAD },
225 { efTPS, NULL, NULL, ffREAD },
226 { efNDX, NULL, NULL, ffOPTRD },
227 { efXVG, NULL, "gyrate", ffWRITE },
228 { efXVG, "-acf", "moi-acf", ffOPTWR },
230 #define NFILE asize(fnm)
235 ppa = add_acf_pargs(&npargs, pa);
237 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW | PCA_BE_NICE,
238 NFILE, fnm, npargs, ppa, asize(desc), desc, 0, NULL, &oenv))
242 bACF = opt2bSet("-acf", NFILE, fnm);
243 if (bACF && nmol != 1)
245 gmx_fatal(FARGS, "Can only do acf with nmol=1");
247 bRot = bRot || bMOI || bACF;
254 printf("Will rotate system along principal axes\n");
255 snew(moi_trans, DIM);
259 printf("Will print moments of inertia\n");
264 printf("Will print radius normalised by charge\n");
267 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &x, NULL, box, TRUE);
268 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
270 if (nmol > gnx || gnx % nmol != 0)
272 gmx_fatal(FARGS, "The number of atoms in the group (%d) is not a multiple of nmol (%d)", gnx, nmol);
276 natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
283 out = xvgropen(ftp2fn(efXVG, NFILE, fnm),
284 "Radius of Charge", "Time (ps)", "Rg (nm)", oenv);
288 out = xvgropen(ftp2fn(efXVG, NFILE, fnm),
289 "Moments of inertia", "Time (ps)", "I (a.m.u. nm\\S2\\N)", oenv);
293 out = xvgropen(ftp2fn(efXVG, NFILE, fnm),
294 "Radius of gyration", "Time (ps)", "Rg (nm)", oenv);
298 xvgr_legend(out, NLEG, legI, oenv);
304 if (output_env_get_print_xvgr_codes(oenv))
306 fprintf(out, "@ subtitle \"Axes are principal component axes\"\n");
309 xvgr_legend(out, NLEG, leg, oenv);
313 gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms);
319 gmx_rmpbc_copy(gpbc, natoms, box, x, x_s);
324 for (mol = 0; mol < nmol; mol++)
326 tm = sub_xcm(nz == 0 ? x_s : x, nam, index+mol*nam, top.atoms.atom, xcm, bQ);
329 gyro += calc_gyro(x_s, nam, index+mol*nam, top.atoms.atom,
330 tm, gvec1, d1, bQ, bRot, bMOI, trans);
334 calc_gyro_z(x, box, nam, index+mol*nam, top.atoms.atom, nz, t, out);
336 rvec_inc(gvec, gvec1);
342 svmul(1.0/nmol, gvec, gvec);
343 svmul(1.0/nmol, d, d);
352 max_moi += delta_moi;
353 for (m = 0; (m < DIM); m++)
355 srenew(moi_trans[m], max_moi*DIM);
358 for (m = 0; (m < DIM); m++)
360 copy_rvec(trans[m], moi_trans[m]+DIM*j);
362 fprintf(out, "%10g %10g %10g %10g %10g\n",
363 t, gyro, d[XX], d[YY], d[ZZ]);
367 fprintf(out, "%10g %10g %10g %10g %10g\n",
368 t, gyro, gvec[XX], gvec[YY], gvec[ZZ]);
373 while (read_next_x(oenv, status, &t, x, box));
377 gmx_rmpbc_done(gpbc);
384 int mode = eacVector;
386 do_autocorr(opt2fn("-acf", NFILE, fnm), oenv,
387 "Moment of inertia vector ACF",
388 j, 3, moi_trans, (t-t0)/j, mode, FALSE);
389 do_view(oenv, opt2fn("-acf", NFILE, fnm), "-nxy");
392 do_view(oenv, ftp2fn(efXVG, NFILE, fnm), "-nxy");