50ec5bee3b111b20d7fb3107e45280fdb319e318
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_gyrate.cpp
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,2016,2017, 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 <cmath>
40 #include <cstring>
41
42 #include "gromacs/commandline/pargs.h"
43 #include "gromacs/commandline/viewit.h"
44 #include "gromacs/correlationfunctions/autocorr.h"
45 #include "gromacs/fileio/confio.h"
46 #include "gromacs/fileio/trxio.h"
47 #include "gromacs/fileio/xvgr.h"
48 #include "gromacs/gmxana/gmx_ana.h"
49 #include "gromacs/gmxana/gstat.h"
50 #include "gromacs/gmxana/princ.h"
51 #include "gromacs/math/functions.h"
52 #include "gromacs/math/utilities.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/pbcutil/rmpbc.h"
55 #include "gromacs/topology/index.h"
56 #include "gromacs/topology/topology.h"
57 #include "gromacs/utility/arraysize.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/futil.h"
60 #include "gromacs/utility/smalloc.h"
61
62 static real calc_gyro(rvec x[], int gnx, int 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] = std::sqrt(d[m]/tm);
80         }
81         /* rotate_atoms(gnx,index,x,trans); */
82     }
83     clear_rvec(comp);
84     for (i = 0; (i < gnx); i++)
85     {
86         ii = index[i];
87         if (bQ)
88         {
89             m0 = std::abs(atom[ii].q);
90         }
91         else
92         {
93             m0 = atom[ii].m;
94         }
95         for (m = 0; (m < DIM); m++)
96         {
97             dx2      = x[ii][m]*x[ii][m];
98             comp[m] += dx2*m0;
99         }
100     }
101     gyro = comp[XX]+comp[YY]+comp[ZZ];
102
103     for (m = 0; (m < DIM); m++)
104     {
105         gvec[m] = std::sqrt((gyro-comp[m])/tm);
106     }
107
108     return std::sqrt(gyro/tm);
109 }
110
111 static void calc_gyro_z(rvec x[], matrix box,
112                         int gnx, int index[], t_atom atom[],
113                         int nz, real time, FILE *out)
114 {
115     static dvec   *inertia = nullptr;
116     static double *tm      = nullptr;
117     int            i, ii, j, zi;
118     real           zf, w, sdet, e1, e2;
119
120     if (inertia == nullptr)
121     {
122         snew(inertia, nz);
123         snew(tm, nz);
124     }
125
126     for (i = 0; i < nz; i++)
127     {
128         clear_dvec(inertia[i]);
129         tm[i] = 0;
130     }
131
132     for (i = 0; (i < gnx); i++)
133     {
134         ii = index[i];
135         zf = nz*x[ii][ZZ]/box[ZZ][ZZ];
136         if (zf >= nz)
137         {
138             zf -= nz;
139         }
140         if (zf < 0)
141         {
142             zf += nz;
143         }
144         for (j = 0; j < 2; j++)
145         {
146             zi = static_cast<int>(zf + j);
147             if (zi == nz)
148             {
149                 zi = 0;
150             }
151             w               = atom[ii].m*(1 + std::cos(M_PI*(zf - zi)));
152             inertia[zi][0] += w*gmx::square(x[ii][YY]);
153             inertia[zi][1] += w*gmx::square(x[ii][XX]);
154             inertia[zi][2] -= w*x[ii][XX]*x[ii][YY];
155             tm[zi]         += w;
156         }
157     }
158     fprintf(out, "%10g", time);
159     for (j = 0; j < nz; j++)
160     {
161         for (i = 0; i < 3; i++)
162         {
163             inertia[j][i] /= tm[j];
164         }
165         sdet = std::sqrt(gmx::square(inertia[j][0] - inertia[j][1]) + 4*gmx::square(inertia[j][2]));
166         e1   = std::sqrt(0.5*(inertia[j][0] + inertia[j][1] + sdet));
167         e2   = std::sqrt(0.5*(inertia[j][0] + inertia[j][1] - sdet));
168         fprintf(out, " %5.3f %5.3f", e1, e2);
169     }
170     fprintf(out, "\n");
171 }
172
173
174 int gmx_gyrate(int argc, char *argv[])
175 {
176     const char       *desc[] = {
177         "[THISMODULE] computes the radius of gyration of a molecule",
178         "and the radii of gyration about the [IT]x[it]-, [IT]y[it]- and [IT]z[it]-axes,",
179         "as a function of time. The atoms are explicitly mass weighted.[PAR]",
180         "The axis components corresponds to the mass-weighted root-mean-square",
181         "of the radii components orthogonal to each axis, for example:[PAR]",
182         "Rg(x) = sqrt((sum_i m_i (R_i(y)^2 + R_i(z)^2))/(sum_i m_i)).[PAR]",
183         "With the [TT]-nmol[tt] option the radius of gyration will be calculated",
184         "for multiple molecules by splitting the analysis group in equally",
185         "sized parts.[PAR]",
186         "With the option [TT]-nz[tt] 2D radii of gyration in the [IT]x-y[it] plane",
187         "of slices along the [IT]z[it]-axis are calculated."
188     };
189     static int        nmol = 1, nz = 0;
190     static gmx_bool   bQ   = FALSE, bRot = FALSE, bMOI = FALSE;
191     t_pargs           pa[] = {
192         { "-nmol", FALSE, etINT, {&nmol},
193           "The number of molecules to analyze" },
194         { "-q", FALSE, etBOOL, {&bQ},
195           "Use absolute value of the charge of an atom as weighting factor instead of mass" },
196         { "-p", FALSE, etBOOL, {&bRot},
197           "Calculate the radii of gyration about the principal axes." },
198         { "-moi", FALSE, etBOOL, {&bMOI},
199           "Calculate the moments of inertia (defined by the principal axes)." },
200         { "-nz", FALSE, etINT, {&nz},
201           "Calculate the 2D radii of gyration of this number of slices along the z-axis" },
202     };
203     FILE             *out;
204     t_trxstatus      *status;
205     t_topology        top;
206     int               ePBC;
207     rvec             *x, *x_s;
208     rvec              xcm, gvec, gvec1;
209     matrix            box, trans;
210     gmx_bool          bACF;
211     real            **moi_trans = nullptr;
212     int               max_moi   = 0, delta_moi = 100;
213     rvec              d, d1; /* eigenvalues of inertia tensor */
214     real              t, t0, tm, gyro;
215     int               natoms;
216     char             *grpname;
217     int               j, m, gnx, nam, mol;
218     int              *index;
219     gmx_output_env_t *oenv;
220     gmx_rmpbc_t       gpbc   = nullptr;
221     const char       *leg[]  = { "Rg", "Rg\\sX\\N", "Rg\\sY\\N", "Rg\\sZ\\N" };
222     const char       *legI[] = { "Itot", "I1", "I2", "I3" };
223 #define NLEG asize(leg)
224     t_filenm          fnm[] = {
225         { efTRX, "-f",   nullptr,       ffREAD },
226         { efTPS, nullptr,   nullptr,       ffREAD },
227         { efNDX, nullptr,   nullptr,       ffOPTRD },
228         { efXVG, nullptr,   "gyrate",   ffWRITE },
229         { efXVG, "-acf", "moi-acf",  ffOPTWR },
230     };
231 #define NFILE asize(fnm)
232     int               npargs;
233     t_pargs          *ppa;
234
235     npargs = asize(pa);
236     ppa    = add_acf_pargs(&npargs, pa);
237
238     if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW,
239                            NFILE, fnm, npargs, ppa, asize(desc), desc, 0, nullptr, &oenv))
240     {
241         sfree(ppa);
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), &top, &ePBC, &x, nullptr, 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 (total and around axes)", "Time (ps)", "Rg (nm)", oenv);
287     }
288     else if (bMOI)
289     {
290         out = xvgropen(ftp2fn(efXVG, NFILE, fnm),
291                        "Moments of inertia (total and around axes)", "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 (total and around axes)", "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(gvec1);
326         clear_rvec(d);
327         clear_rvec(d1);
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_trx(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 }