Further copyrite.h cleanup.
[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 "futil.h"
50 #include "statutil.h"
51 #include "index.h"
52 #include "mshift.h"
53 #include "xvgr.h"
54 #include "princ.h"
55 #include "rmpbc.h"
56 #include "txtdump.h"
57 #include "tpxio.h"
58 #include "gstat.h"
59 #include "gmx_ana.h"
60
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 int gmx_gyrate(int argc, char *argv[])
177 {
178     const char     *desc[] = {
179         "[TT]g_gyrate[tt] computes the radius of gyration of a molecule",
180         "and the radii of gyration about the [IT]x[it]-, [IT]y[it]- and [IT]z[it]-axes,",
181         "as a function of time. The atoms are explicitly mass weighted.[PAR]",
182         "With the [TT]-nmol[tt] option the radius of gyration will be calculated",
183         "for multiple molecules by splitting the analysis group in equally",
184         "sized parts.[PAR]",
185         "With the option [TT]-nz[tt] 2D radii of gyration in the [IT]x-y[it] plane",
186         "of slices along the [IT]z[it]-axis are calculated."
187     };
188     static int      nmol = 1, nz = 0;
189     static gmx_bool bQ   = FALSE, bRot = FALSE, bMOI = FALSE;
190     t_pargs         pa[] = {
191         { "-nmol", FALSE, etINT, {&nmol},
192           "The number of molecules to analyze" },
193         { "-q", FALSE, etBOOL, {&bQ},
194           "Use absolute value of the charge of an atom as weighting factor instead of mass" },
195         { "-p", FALSE, etBOOL, {&bRot},
196           "Calculate the radii of gyration about the principal axes." },
197         { "-moi", FALSE, etBOOL, {&bMOI},
198           "Calculate the moments of inertia (defined by the principal axes)." },
199         { "-nz", FALSE, etINT, {&nz},
200           "Calculate the 2D radii of gyration of this number of slices along the z-axis" },
201     };
202     FILE           *out;
203     t_trxstatus    *status;
204     t_topology      top;
205     int             ePBC;
206     rvec           *x, *x_s;
207     rvec            xcm, gvec, gvec1;
208     matrix          box, trans;
209     gmx_bool        bACF;
210     real          **moi_trans = NULL;
211     int             max_moi   = 0, delta_moi = 100;
212     rvec            d, d1; /* eigenvalues of inertia tensor */
213     real            t, t0, tm, gyro;
214     int             natoms;
215     char           *grpname, title[256];
216     int             i, j, m, gnx, nam, mol;
217     atom_id        *index;
218     output_env_t    oenv;
219     gmx_rmpbc_t     gpbc   = NULL;
220     const char     *leg[]  = { "Rg", "RgX", "RgY", "RgZ" };
221     const char     *legI[] = { "Itot", "I1", "I2", "I3" };
222 #define NLEG asize(leg)
223     t_filenm        fnm[] = {
224         { efTRX, "-f",   NULL,       ffREAD },
225         { efTPS, NULL,   NULL,       ffREAD },
226         { efNDX, NULL,   NULL,       ffOPTRD },
227         { efXVG, NULL,   "gyrate",   ffWRITE },
228         { efXVG, "-acf", "moi-acf",  ffOPTWR },
229     };
230 #define NFILE asize(fnm)
231     int             npargs;
232     t_pargs        *ppa;
233
234     npargs = asize(pa);
235     ppa    = add_acf_pargs(&npargs, pa);
236
237     parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW | PCA_BE_NICE,
238                       NFILE, fnm, npargs, ppa, asize(desc), desc, 0, NULL, &oenv);
239     bACF = opt2bSet("-acf", NFILE, fnm);
240     if (bACF && nmol != 1)
241     {
242         gmx_fatal(FARGS, "Can only do acf with nmol=1");
243     }
244     bRot = bRot || bMOI || bACF;
245     /*
246        if (nz > 0)
247        bMOI = TRUE;
248      */
249     if (bRot)
250     {
251         printf("Will rotate system along principal axes\n");
252         snew(moi_trans, DIM);
253     }
254     if (bMOI)
255     {
256         printf("Will print moments of inertia\n");
257         bQ = FALSE;
258     }
259     if (bQ)
260     {
261         printf("Will print radius normalised by charge\n");
262     }
263
264     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &x, NULL, box, TRUE);
265     get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
266
267     if (nmol > gnx || gnx % nmol != 0)
268     {
269         gmx_fatal(FARGS, "The number of atoms in the group (%d) is not a multiple of nmol (%d)", gnx, nmol);
270     }
271     nam = gnx/nmol;
272
273     natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
274     snew(x_s, natoms);
275
276     j  = 0;
277     t0 = t;
278     if (bQ)
279     {
280         out = xvgropen(ftp2fn(efXVG, NFILE, fnm),
281                        "Radius of Charge", "Time (ps)", "Rg (nm)", oenv);
282     }
283     else if (bMOI)
284     {
285         out = xvgropen(ftp2fn(efXVG, NFILE, fnm),
286                        "Moments of inertia", "Time (ps)", "I (a.m.u. nm\\S2\\N)", oenv);
287     }
288     else
289     {
290         out = xvgropen(ftp2fn(efXVG, NFILE, fnm),
291                        "Radius of gyration", "Time (ps)", "Rg (nm)", oenv);
292     }
293     if (bMOI)
294     {
295         xvgr_legend(out, NLEG, legI, oenv);
296     }
297     else
298     {
299         if (bRot)
300         {
301             if (output_env_get_print_xvgr_codes(oenv))
302             {
303                 fprintf(out, "@ subtitle \"Axes are principal component axes\"\n");
304             }
305         }
306         xvgr_legend(out, NLEG, leg, oenv);
307     }
308     if (nz == 0)
309     {
310         gpbc = gmx_rmpbc_init(&top.idef, ePBC, natoms, box);
311     }
312     do
313     {
314         if (nz == 0)
315         {
316             gmx_rmpbc_copy(gpbc, natoms, box, x, x_s);
317         }
318         gyro = 0;
319         clear_rvec(gvec);
320         clear_rvec(d);
321         for (mol = 0; mol < nmol; mol++)
322         {
323             tm    = sub_xcm(nz == 0 ? x_s : x, nam, index+mol*nam, top.atoms.atom, xcm, bQ);
324             if (nz == 0)
325             {
326                 gyro += calc_gyro(x_s, nam, index+mol*nam, top.atoms.atom,
327                                   tm, gvec1, d1, bQ, bRot, bMOI, trans);
328             }
329             else
330             {
331                 calc_gyro_z(x, box, nam, index+mol*nam, top.atoms.atom, nz, t, out);
332             }
333             rvec_inc(gvec, gvec1);
334             rvec_inc(d, d1);
335         }
336         if (nmol > 0)
337         {
338             gyro /= nmol;
339             svmul(1.0/nmol, gvec, gvec);
340             svmul(1.0/nmol, d, d);
341         }
342
343         if (nz == 0)
344         {
345             if (bRot)
346             {
347                 if (j >= max_moi)
348                 {
349                     max_moi += delta_moi;
350                     for (m = 0; (m < DIM); m++)
351                     {
352                         srenew(moi_trans[m], max_moi*DIM);
353                     }
354                 }
355                 for (m = 0; (m < DIM); m++)
356                 {
357                     copy_rvec(trans[m], moi_trans[m]+DIM*j);
358                 }
359                 fprintf(out, "%10g  %10g  %10g  %10g  %10g\n",
360                         t, gyro, d[XX], d[YY], d[ZZ]);
361             }
362             else
363             {
364                 fprintf(out, "%10g  %10g  %10g  %10g  %10g\n",
365                         t, gyro, gvec[XX], gvec[YY], gvec[ZZ]);
366             }
367         }
368         j++;
369     }
370     while (read_next_x(oenv, status, &t, natoms, x, box));
371     close_trj(status);
372     if (nz == 0)
373     {
374         gmx_rmpbc_done(gpbc);
375     }
376
377     ffclose(out);
378
379     if (bACF)
380     {
381         int mode = eacVector;
382
383         do_autocorr(opt2fn("-acf", NFILE, fnm), oenv,
384                     "Moment of inertia vector ACF",
385                     j, 3, moi_trans, (t-t0)/j, mode, FALSE);
386         do_view(oenv, opt2fn("-acf", NFILE, fnm), "-nxy");
387     }
388
389     do_view(oenv, ftp2fn(efXVG, NFILE, fnm), "-nxy");
390
391     return 0;
392 }