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, 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.
42 #include "gromacs/commandline/pargs.h"
43 #include "gromacs/fileio/tpxio.h"
44 #include "gromacs/fileio/trxio.h"
45 #include "gromacs/fileio/xvgr.h"
46 #include "gromacs/gmxana/gmx_ana.h"
47 #include "gromacs/gmxana/gstat.h"
48 #include "gromacs/gmxana/princ.h"
49 #include "gromacs/legacyheaders/macros.h"
50 #include "gromacs/legacyheaders/txtdump.h"
51 #include "gromacs/legacyheaders/typedefs.h"
52 #include "gromacs/legacyheaders/viewit.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/pbcutil/rmpbc.h"
55 #include "gromacs/topology/index.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/futil.h"
58 #include "gromacs/utility/smalloc.h"
60 real calc_gyro(rvec x[], int gnx, atom_id index[], t_atom atom[], real tm,
61 rvec gvec, rvec d, gmx_bool bQ, gmx_bool bRot, gmx_bool bMOI, matrix trans)
64 real gyro, dx2, m0, Itot;
69 principal_comp(gnx, index, atom, x, trans, d);
75 for (m = 0; (m < DIM); m++)
80 pr_rvecs(stderr, 0, "trans", trans, DIM);
82 /* rotate_atoms(gnx,index,x,trans); */
85 for (i = 0; (i < gnx); i++)
90 m0 = fabs(atom[ii].q);
96 for (m = 0; (m < DIM); m++)
98 dx2 = x[ii][m]*x[ii][m];
102 gyro = comp[XX]+comp[YY]+comp[ZZ];
104 for (m = 0; (m < DIM); m++)
106 gvec[m] = sqrt((gyro-comp[m])/tm);
109 return sqrt(gyro/tm);
112 void calc_gyro_z(rvec x[], matrix box,
113 int gnx, atom_id index[], t_atom atom[],
114 int nz, real time, FILE *out)
116 static dvec *inertia = NULL;
117 static double *tm = NULL;
119 real zf, w, sdet, e1, e2;
127 for (i = 0; i < nz; i++)
129 clear_dvec(inertia[i]);
133 for (i = 0; (i < gnx); i++)
136 zf = nz*x[ii][ZZ]/box[ZZ][ZZ];
145 for (j = 0; j < 2; j++)
152 w = atom[ii].m*(1 + cos(M_PI*(zf - zi)));
153 inertia[zi][0] += w*sqr(x[ii][YY]);
154 inertia[zi][1] += w*sqr(x[ii][XX]);
155 inertia[zi][2] -= w*x[ii][XX]*x[ii][YY];
159 fprintf(out, "%10g", time);
160 for (j = 0; j < nz; j++)
162 for (i = 0; i < 3; i++)
164 inertia[j][i] /= tm[j];
166 sdet = sqrt(sqr(inertia[j][0] - inertia[j][1]) + 4*sqr(inertia[j][2]));
167 e1 = sqrt(0.5*(inertia[j][0] + inertia[j][1] + sdet));
168 e2 = sqrt(0.5*(inertia[j][0] + inertia[j][1] - sdet));
169 fprintf(out, " %5.3f %5.3f", e1, e2);
174 int gmx_gyrate(int argc, char *argv[])
176 const char *desc[] = {
177 "[THISMODULE] computes the radius of gyration of a molecule",
178 "and the radii of gyration about the [IT]x[it]-, [IT]y[it]- and [IT]z[it]-axes,",
179 "as a function of time. The atoms are explicitly mass weighted.[PAR]",
180 "With the [TT]-nmol[tt] option the radius of gyration will be calculated",
181 "for multiple molecules by splitting the analysis group in equally",
183 "With the option [TT]-nz[tt] 2D radii of gyration in the [IT]x-y[it] plane",
184 "of slices along the [IT]z[it]-axis are calculated."
186 static int nmol = 1, nz = 0;
187 static gmx_bool bQ = FALSE, bRot = FALSE, bMOI = FALSE;
189 { "-nmol", FALSE, etINT, {&nmol},
190 "The number of molecules to analyze" },
191 { "-q", FALSE, etBOOL, {&bQ},
192 "Use absolute value of the charge of an atom as weighting factor instead of mass" },
193 { "-p", FALSE, etBOOL, {&bRot},
194 "Calculate the radii of gyration about the principal axes." },
195 { "-moi", FALSE, etBOOL, {&bMOI},
196 "Calculate the moments of inertia (defined by the principal axes)." },
197 { "-nz", FALSE, etINT, {&nz},
198 "Calculate the 2D radii of gyration of this number of slices along the z-axis" },
205 rvec xcm, gvec, gvec1;
208 real **moi_trans = NULL;
209 int max_moi = 0, delta_moi = 100;
210 rvec d, d1; /* eigenvalues of inertia tensor */
211 real t, t0, tm, gyro;
213 char *grpname, title[256];
214 int i, j, m, gnx, nam, mol;
217 gmx_rmpbc_t gpbc = NULL;
218 const char *leg[] = { "Rg", "RgX", "RgY", "RgZ" };
219 const char *legI[] = { "Itot", "I1", "I2", "I3" };
220 #define NLEG asize(leg)
222 { efTRX, "-f", NULL, ffREAD },
223 { efTPS, NULL, NULL, ffREAD },
224 { efNDX, NULL, NULL, ffOPTRD },
225 { efXVG, NULL, "gyrate", ffWRITE },
226 { efXVG, "-acf", "moi-acf", ffOPTWR },
228 #define NFILE asize(fnm)
233 ppa = add_acf_pargs(&npargs, pa);
235 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW,
236 NFILE, fnm, npargs, ppa, asize(desc), desc, 0, NULL, &oenv))
240 bACF = opt2bSet("-acf", NFILE, fnm);
241 if (bACF && nmol != 1)
243 gmx_fatal(FARGS, "Can only do acf with nmol=1");
245 bRot = bRot || bMOI || bACF;
252 printf("Will rotate system along principal axes\n");
253 snew(moi_trans, DIM);
257 printf("Will print moments of inertia\n");
262 printf("Will print radius normalised by charge\n");
265 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &x, NULL, box, TRUE);
266 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
268 if (nmol > gnx || gnx % nmol != 0)
270 gmx_fatal(FARGS, "The number of atoms in the group (%d) is not a multiple of nmol (%d)", gnx, nmol);
274 natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
281 out = xvgropen(ftp2fn(efXVG, NFILE, fnm),
282 "Radius of Charge", "Time (ps)", "Rg (nm)", oenv);
286 out = xvgropen(ftp2fn(efXVG, NFILE, fnm),
287 "Moments of inertia", "Time (ps)", "I (a.m.u. nm\\S2\\N)", oenv);
291 out = xvgropen(ftp2fn(efXVG, NFILE, fnm),
292 "Radius of gyration", "Time (ps)", "Rg (nm)", oenv);
296 xvgr_legend(out, NLEG, legI, oenv);
302 if (output_env_get_print_xvgr_codes(oenv))
304 fprintf(out, "@ subtitle \"Axes are principal component axes\"\n");
307 xvgr_legend(out, NLEG, leg, oenv);
311 gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms);
317 gmx_rmpbc_copy(gpbc, natoms, box, x, x_s);
322 for (mol = 0; mol < nmol; mol++)
324 tm = sub_xcm(nz == 0 ? x_s : x, nam, index+mol*nam, top.atoms.atom, xcm, bQ);
327 gyro += calc_gyro(x_s, nam, index+mol*nam, top.atoms.atom,
328 tm, gvec1, d1, bQ, bRot, bMOI, trans);
332 calc_gyro_z(x, box, nam, index+mol*nam, top.atoms.atom, nz, t, out);
334 rvec_inc(gvec, gvec1);
340 svmul(1.0/nmol, gvec, gvec);
341 svmul(1.0/nmol, d, d);
350 max_moi += delta_moi;
351 for (m = 0; (m < DIM); m++)
353 srenew(moi_trans[m], max_moi*DIM);
356 for (m = 0; (m < DIM); m++)
358 copy_rvec(trans[m], moi_trans[m]+DIM*j);
360 fprintf(out, "%10g %10g %10g %10g %10g\n",
361 t, gyro, d[XX], d[YY], d[ZZ]);
365 fprintf(out, "%10g %10g %10g %10g %10g\n",
366 t, gyro, gvec[XX], gvec[YY], gvec[ZZ]);
371 while (read_next_x(oenv, status, &t, x, box));
375 gmx_rmpbc_done(gpbc);
382 int mode = eacVector;
384 do_autocorr(opt2fn("-acf", NFILE, fnm), oenv,
385 "Moment of inertia vector ACF",
386 j, 3, moi_trans, (t-t0)/j, mode, FALSE);
387 do_view(oenv, opt2fn("-acf", NFILE, fnm), "-nxy");
390 do_view(oenv, ftp2fn(efXVG, NFILE, fnm), "-nxy");