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