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