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