Merge branch 'release-4-6'
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_gyrate.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
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.
14
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.
19  *
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.
26  *
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.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Green Red Orange Magenta Azure Cyan Skyblue
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <math.h>
40 #include <string.h>
41
42 #include "statutil.h"
43 #include "sysstuff.h"
44 #include "typedefs.h"
45 #include "smalloc.h"
46 #include "macros.h"
47 #include "vec.h"
48 #include "pbc.h"
49 #include "copyrite.h"
50 #include "futil.h"
51 #include "statutil.h"
52 #include "index.h"
53 #include "mshift.h"
54 #include "xvgr.h"
55 #include "princ.h"
56 #include "rmpbc.h"
57 #include "txtdump.h"
58 #include "tpxio.h"
59 #include "gstat.h"
60 #include "gmx_ana.h"
61
62
63 real calc_gyro(rvec x[], int gnx, atom_id index[], t_atom atom[], real tm,
64                rvec gvec, rvec d, gmx_bool bQ, gmx_bool bRot, gmx_bool bMOI, matrix trans)
65 {
66     int    i, ii, m;
67     real   gyro, dx2, m0, Itot;
68     rvec   comp;
69
70     if (bRot)
71     {
72         principal_comp(gnx, index, atom, x, trans, d);
73         Itot = norm(d);
74         if (bMOI)
75         {
76             return Itot;
77         }
78         for (m = 0; (m < DIM); m++)
79         {
80             d[m] = sqrt(d[m]/tm);
81         }
82 #ifdef DEBUG
83         pr_rvecs(stderr, 0, "trans", trans, DIM);
84 #endif
85         /* rotate_atoms(gnx,index,x,trans); */
86     }
87     clear_rvec(comp);
88     for (i = 0; (i < gnx); i++)
89     {
90         ii = index[i];
91         if (bQ)
92         {
93             m0 = fabs(atom[ii].q);
94         }
95         else
96         {
97             m0 = atom[ii].m;
98         }
99         for (m = 0; (m < DIM); m++)
100         {
101             dx2      = x[ii][m]*x[ii][m];
102             comp[m] += dx2*m0;
103         }
104     }
105     gyro = comp[XX]+comp[YY]+comp[ZZ];
106
107     for (m = 0; (m < DIM); m++)
108     {
109         gvec[m] = sqrt((gyro-comp[m])/tm);
110     }
111
112     return sqrt(gyro/tm);
113 }
114
115 void calc_gyro_z(rvec x[], matrix box,
116                  int gnx, atom_id index[], t_atom atom[],
117                  int nz, real time, FILE *out)
118 {
119     static dvec   *inertia = NULL;
120     static double *tm      = NULL;
121     int            i, ii, j, zi;
122     real           zf, w, sdet, e1, e2;
123
124     if (inertia == NULL)
125     {
126         snew(inertia, nz);
127         snew(tm, nz);
128     }
129
130     for (i = 0; i < nz; i++)
131     {
132         clear_dvec(inertia[i]);
133         tm[i] = 0;
134     }
135
136     for (i = 0; (i < gnx); i++)
137     {
138         ii = index[i];
139         zf = nz*x[ii][ZZ]/box[ZZ][ZZ];
140         if (zf >= nz)
141         {
142             zf -= nz;
143         }
144         if (zf < 0)
145         {
146             zf += nz;
147         }
148         for (j = 0; j < 2; j++)
149         {
150             zi = zf + j;
151             if (zi == nz)
152             {
153                 zi = 0;
154             }
155             w               = atom[ii].m*(1 + cos(M_PI*(zf - zi)));
156             inertia[zi][0] += w*sqr(x[ii][YY]);
157             inertia[zi][1] += w*sqr(x[ii][XX]);
158             inertia[zi][2] -= w*x[ii][XX]*x[ii][YY];
159             tm[zi]         += w;
160         }
161     }
162     fprintf(out, "%10g", time);
163     for (j = 0; j < nz; j++)
164     {
165         for (i = 0; i < 3; i++)
166         {
167             inertia[j][i] /= tm[j];
168         }
169         sdet = sqrt(sqr(inertia[j][0] - inertia[j][1]) + 4*sqr(inertia[j][2]));
170         e1   = sqrt(0.5*(inertia[j][0] + inertia[j][1] + sdet));
171         e2   = sqrt(0.5*(inertia[j][0] + inertia[j][1] - sdet));
172         fprintf(out, " %5.3f %5.3f", e1, e2);
173     }
174     fprintf(out, "\n");
175 }
176
177 int gmx_gyrate(int argc, char *argv[])
178 {
179     const char     *desc[] = {
180         "[TT]g_gyrate[tt] computes the radius of gyration of a group of atoms",
181         "and the radii of gyration about the [IT]x[it]-, [IT]y[it]- and [IT]z[it]-axes,",
182         "as a function of time. The atoms are explicitly mass weighted.[PAR]",
183         "With the [TT]-nmol[tt] option the radius of gyration will be calculated",
184         "for multiple molecules by splitting the analysis group in equally",
185         "sized parts.[PAR]",
186         "With the option [TT]-nz[tt] 2D radii of gyration in the [IT]x-y[it] plane",
187         "of slices along the [IT]z[it]-axis are calculated."
188     };
189     static int      nmol = 1, nz = 0;
190     static gmx_bool bQ   = FALSE, bRot = FALSE, bMOI = FALSE;
191     t_pargs         pa[] = {
192         { "-nmol", FALSE, etINT, {&nmol},
193           "The number of molecules to analyze" },
194         { "-q", FALSE, etBOOL, {&bQ},
195           "Use absolute value of the charge of an atom as weighting factor instead of mass" },
196         { "-p", FALSE, etBOOL, {&bRot},
197           "Calculate the radii of gyration about the principal axes." },
198         { "-moi", FALSE, etBOOL, {&bMOI},
199           "Calculate the moments of inertia (defined by the principal axes)." },
200         { "-nz", FALSE, etINT, {&nz},
201           "Calculate the 2D radii of gyration of this number of slices along the z-axis" },
202     };
203     FILE           *out;
204     t_trxstatus    *status;
205     t_topology      top;
206     int             ePBC;
207     rvec           *x, *x_s;
208     rvec            xcm, gvec, gvec1;
209     matrix          box, trans;
210     gmx_bool        bACF;
211     real          **moi_trans = NULL;
212     int             max_moi   = 0, delta_moi = 100;
213     rvec            d, d1; /* eigenvalues of inertia tensor */
214     real            t, t0, tm, gyro;
215     int             natoms;
216     char           *grpname, title[256];
217     int             i, j, m, gnx, nam, mol;
218     atom_id        *index;
219     output_env_t    oenv;
220     gmx_rmpbc_t     gpbc   = NULL;
221     const char     *leg[]  = { "Rg", "RgX", "RgY", "RgZ" };
222     const char     *legI[] = { "Itot", "I1", "I2", "I3" };
223 #define NLEG asize(leg)
224     t_filenm        fnm[] = {
225         { efTRX, "-f",   NULL,       ffREAD },
226         { efTPS, NULL,   NULL,       ffREAD },
227         { efNDX, NULL,   NULL,       ffOPTRD },
228         { efXVG, NULL,   "gyrate",   ffWRITE },
229         { efXVG, "-acf", "moi-acf",  ffOPTWR },
230     };
231 #define NFILE asize(fnm)
232     int             npargs;
233     t_pargs        *ppa;
234
235     npargs = asize(pa);
236     ppa    = add_acf_pargs(&npargs, pa);
237
238     parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW | PCA_BE_NICE,
239                       NFILE, fnm, npargs, ppa, asize(desc), desc, 0, NULL, &oenv);
240     bACF = opt2bSet("-acf", NFILE, fnm);
241     if (bACF && nmol != 1)
242     {
243         gmx_fatal(FARGS, "Can only do acf with nmol=1");
244     }
245     bRot = bRot || bMOI || bACF;
246     /*
247        if (nz > 0)
248        bMOI = TRUE;
249      */
250     if (bRot)
251     {
252         printf("Will rotate system along principal axes\n");
253         snew(moi_trans, DIM);
254     }
255     if (bMOI)
256     {
257         printf("Will print moments of inertia\n");
258         bQ = FALSE;
259     }
260     if (bQ)
261     {
262         printf("Will print radius normalised by charge\n");
263     }
264
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);
267
268     if (nmol > gnx || gnx % nmol != 0)
269     {
270         gmx_fatal(FARGS, "The number of atoms in the group (%d) is not a multiple of nmol (%d)", gnx, nmol);
271     }
272     nam = gnx/nmol;
273
274     natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
275     snew(x_s, natoms);
276
277     j  = 0;
278     t0 = t;
279     if (bQ)
280     {
281         out = xvgropen(ftp2fn(efXVG, NFILE, fnm),
282                        "Radius of Charge", "Time (ps)", "Rg (nm)", oenv);
283     }
284     else if (bMOI)
285     {
286         out = xvgropen(ftp2fn(efXVG, NFILE, fnm),
287                        "Moments of inertia", "Time (ps)", "I (a.m.u. nm\\S2\\N)", oenv);
288     }
289     else
290     {
291         out = xvgropen(ftp2fn(efXVG, NFILE, fnm),
292                        "Radius of gyration", "Time (ps)", "Rg (nm)", oenv);
293     }
294     if (bMOI)
295     {
296         xvgr_legend(out, NLEG, legI, oenv);
297     }
298     else
299     {
300         if (bRot)
301         {
302             if (output_env_get_print_xvgr_codes(oenv))
303             {
304                 fprintf(out, "@ subtitle \"Axes are principal component axes\"\n");
305             }
306         }
307         xvgr_legend(out, NLEG, leg, oenv);
308     }
309     if (nz == 0)
310     {
311         gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms, box);
312     }
313     do
314     {
315         if (nz == 0)
316         {
317             gmx_rmpbc_copy(gpbc, natoms, box, x, x_s);
318         }
319         gyro = 0;
320         clear_rvec(gvec);
321         clear_rvec(d);
322         for (mol = 0; mol < nmol; mol++)
323         {
324             tm    = sub_xcm(nz == 0 ? x_s : x, nam, index+mol*nam, top.atoms.atom, xcm, bQ);
325             if (nz == 0)
326             {
327                 gyro += calc_gyro(x_s, nam, index+mol*nam, top.atoms.atom,
328                                   tm, gvec1, d1, bQ, bRot, bMOI, trans);
329             }
330             else
331             {
332                 calc_gyro_z(x, box, nam, index+mol*nam, top.atoms.atom, nz, t, out);
333             }
334             rvec_inc(gvec, gvec1);
335             rvec_inc(d, d1);
336         }
337         if (nmol > 0)
338         {
339             gyro /= nmol;
340             svmul(1.0/nmol, gvec, gvec);
341             svmul(1.0/nmol, d, d);
342         }
343
344         if (nz == 0)
345         {
346             if (bRot)
347             {
348                 if (j >= max_moi)
349                 {
350                     max_moi += delta_moi;
351                     for (m = 0; (m < DIM); m++)
352                     {
353                         srenew(moi_trans[m], max_moi*DIM);
354                     }
355                 }
356                 for (m = 0; (m < DIM); m++)
357                 {
358                     copy_rvec(trans[m], moi_trans[m]+DIM*j);
359                 }
360                 fprintf(out, "%10g  %10g  %10g  %10g  %10g\n",
361                         t, gyro, d[XX], d[YY], d[ZZ]);
362             }
363             else
364             {
365                 fprintf(out, "%10g  %10g  %10g  %10g  %10g\n",
366                         t, gyro, gvec[XX], gvec[YY], gvec[ZZ]);
367             }
368         }
369         j++;
370     }
371     while (read_next_x(oenv, status, &t, natoms, x, box));
372     close_trj(status);
373     if (nz == 0)
374     {
375         gmx_rmpbc_done(gpbc);
376     }
377
378     ffclose(out);
379
380     if (bACF)
381     {
382         int mode = eacVector;
383
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");
388     }
389
390     do_view(oenv, ftp2fn(efXVG, NFILE, fnm), "-nxy");
391
392     thanx(stderr);
393
394     return 0;
395 }