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