Move M_PI definition to math/units.h
[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.
7  * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #include "gmxpre.h"
39
40 #include <cmath>
41 #include <cstring>
42
43 #include "gromacs/commandline/pargs.h"
44 #include "gromacs/commandline/viewit.h"
45 #include "gromacs/correlationfunctions/autocorr.h"
46 #include "gromacs/fileio/confio.h"
47 #include "gromacs/fileio/trxio.h"
48 #include "gromacs/fileio/xvgr.h"
49 #include "gromacs/gmxana/gmx_ana.h"
50 #include "gromacs/gmxana/gstat.h"
51 #include "gromacs/gmxana/princ.h"
52 #include "gromacs/math/functions.h"
53 #include "gromacs/math/units.h"
54 #include "gromacs/math/utilities.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/pbcutil/rmpbc.h"
57 #include "gromacs/topology/index.h"
58 #include "gromacs/topology/topology.h"
59 #include "gromacs/utility/arraysize.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/futil.h"
62 #include "gromacs/utility/smalloc.h"
63
64 static real calc_gyro(rvec     x[],
65                       int      gnx,
66                       int      index[],
67                       t_atom   atom[],
68                       real     tm,
69                       rvec     gvec,
70                       rvec     d,
71                       gmx_bool bQ,
72                       gmx_bool bRot,
73                       gmx_bool bMOI,
74                       matrix   trans)
75 {
76     int  i, ii, m;
77     real gyro, dx2, m0, Itot;
78     rvec comp;
79
80     if (bRot)
81     {
82         principal_comp(gnx, index, atom, x, trans, d);
83         Itot = norm(d);
84         if (bMOI)
85         {
86             return Itot;
87         }
88         for (m = 0; (m < DIM); m++)
89         {
90             d[m] = std::sqrt(d[m] / tm);
91         }
92         /* rotate_atoms(gnx,index,x,trans); */
93     }
94     clear_rvec(comp);
95     for (i = 0; (i < gnx); i++)
96     {
97         ii = index[i];
98         if (bQ)
99         {
100             m0 = std::abs(atom[ii].q);
101         }
102         else
103         {
104             m0 = atom[ii].m;
105         }
106         for (m = 0; (m < DIM); m++)
107         {
108             dx2 = x[ii][m] * x[ii][m];
109             comp[m] += dx2 * m0;
110         }
111     }
112     gyro = comp[XX] + comp[YY] + comp[ZZ];
113
114     for (m = 0; (m < DIM); m++)
115     {
116         gvec[m] = std::sqrt((gyro - comp[m]) / tm);
117     }
118
119     return std::sqrt(gyro / tm);
120 }
121
122 static void calc_gyro_z(rvec x[], matrix box, int gnx, const int index[], t_atom atom[], int nz, real time, FILE* out)
123 {
124     static dvec*   inertia = nullptr;
125     static double* tm      = nullptr;
126     int            i, ii, j, zi;
127     real           zf, w, sdet, e1, e2;
128
129     if (inertia == nullptr)
130     {
131         snew(inertia, nz);
132         snew(tm, nz);
133     }
134
135     for (i = 0; i < nz; i++)
136     {
137         clear_dvec(inertia[i]);
138         tm[i] = 0;
139     }
140
141     for (i = 0; (i < gnx); i++)
142     {
143         ii = index[i];
144         zf = nz * x[ii][ZZ] / box[ZZ][ZZ];
145         if (zf >= nz)
146         {
147             zf -= nz;
148         }
149         if (zf < 0)
150         {
151             zf += nz;
152         }
153         for (j = 0; j < 2; j++)
154         {
155             zi = static_cast<int>(zf + j);
156             if (zi == nz)
157             {
158                 zi = 0;
159             }
160             w = atom[ii].m * (1 + std::cos(M_PI * (zf - zi)));
161             inertia[zi][0] += w * gmx::square(x[ii][YY]);
162             inertia[zi][1] += w * gmx::square(x[ii][XX]);
163             inertia[zi][2] -= w * x[ii][XX] * x[ii][YY];
164             tm[zi] += w;
165         }
166     }
167     fprintf(out, "%10g", time);
168     for (j = 0; j < nz; j++)
169     {
170         for (i = 0; i < 3; i++)
171         {
172             inertia[j][i] /= tm[j];
173         }
174         sdet = std::sqrt(gmx::square(inertia[j][0] - inertia[j][1]) + 4 * gmx::square(inertia[j][2]));
175         e1   = std::sqrt(0.5 * (inertia[j][0] + inertia[j][1] + sdet));
176         e2   = std::sqrt(0.5 * (inertia[j][0] + inertia[j][1] - sdet));
177         fprintf(out, " %5.3f %5.3f", e1, e2);
178     }
179     fprintf(out, "\n");
180 }
181
182
183 int gmx_gyrate(int argc, char* argv[])
184 {
185     const char* desc[] = {
186         "[THISMODULE] computes the radius of gyration of a molecule",
187         "and the radii of gyration about the [IT]x[it]-, [IT]y[it]- and [IT]z[it]-axes,",
188         "as a function of time. The atoms are explicitly mass weighted.[PAR]",
189         "The axis components corresponds to the mass-weighted root-mean-square",
190         "of the radii components orthogonal to each axis, for example:[PAR]",
191         "Rg(x) = sqrt((sum_i m_i (R_i(y)^2 + R_i(z)^2))/(sum_i m_i)).[PAR]",
192         "With the [TT]-nmol[tt] option the radius of gyration will be calculated",
193         "for multiple molecules by splitting the analysis group in equally",
194         "sized parts.[PAR]",
195         "With the option [TT]-nz[tt] 2D radii of gyration in the [IT]x-y[it] plane",
196         "of slices along the [IT]z[it]-axis are calculated."
197     };
198     static int      nmol = 1, nz = 0;
199     static gmx_bool bQ = FALSE, bRot = FALSE, bMOI = FALSE;
200     t_pargs         pa[] = {
201         { "-nmol", FALSE, etINT, { &nmol }, "The number of molecules to analyze" },
202         { "-q",
203           FALSE,
204           etBOOL,
205           { &bQ },
206           "Use absolute value of the charge of an atom as weighting factor instead of mass" },
207         { "-p",
208           FALSE,
209           etBOOL,
210           { &bRot },
211           "Calculate the radii of gyration about the principal axes." },
212         { "-moi",
213           FALSE,
214           etBOOL,
215           { &bMOI },
216           "Calculate the moments of inertia (defined by the principal axes)." },
217         { "-nz",
218           FALSE,
219           etINT,
220           { &nz },
221           "Calculate the 2D radii of gyration of this number of slices along the z-axis" },
222     };
223     FILE*             out;
224     t_trxstatus*      status;
225     t_topology        top;
226     PbcType           pbcType;
227     rvec *            x, *x_s;
228     rvec              xcm, gvec, gvec1;
229     matrix            box, trans;
230     gmx_bool          bACF;
231     real**            moi_trans = nullptr;
232     int               max_moi = 0, delta_moi = 100;
233     rvec              d, d1; /* eigenvalues of inertia tensor */
234     real              t, t0, tm, gyro;
235     int               natoms;
236     char*             grpname;
237     int               j, m, gnx, nam, mol;
238     int*              index;
239     gmx_output_env_t* oenv;
240     gmx_rmpbc_t       gpbc   = nullptr;
241     const char*       leg[]  = { "Rg", "Rg\\sX\\N", "Rg\\sY\\N", "Rg\\sZ\\N" };
242     const char*       legI[] = { "Itot", "I1", "I2", "I3" };
243 #define NLEG asize(leg)
244     t_filenm fnm[] = {
245         { efTRX, "-f", nullptr, ffREAD },      { efTPS, nullptr, nullptr, ffREAD },
246         { efNDX, nullptr, nullptr, ffOPTRD },  { efXVG, nullptr, "gyrate", ffWRITE },
247         { efXVG, "-acf", "moi-acf", ffOPTWR },
248     };
249 #define NFILE asize(fnm)
250     int      npargs;
251     t_pargs* ppa;
252
253     npargs = asize(pa);
254     ppa    = add_acf_pargs(&npargs, pa);
255
256     if (!parse_common_args(
257                 &argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW, NFILE, fnm, npargs, ppa, asize(desc), desc, 0, nullptr, &oenv))
258     {
259         sfree(ppa);
260         return 0;
261     }
262     bACF = opt2bSet("-acf", NFILE, fnm);
263     if (bACF && nmol != 1)
264     {
265         gmx_fatal(FARGS, "Can only do acf with nmol=1");
266     }
267     bRot = bRot || bMOI || bACF;
268     /*
269        if (nz > 0)
270        bMOI = TRUE;
271      */
272     if (bRot)
273     {
274         printf("Will rotate system along principal axes\n");
275         snew(moi_trans, DIM);
276     }
277     if (bMOI)
278     {
279         printf("Will print moments of inertia\n");
280         bQ = FALSE;
281     }
282     if (bQ)
283     {
284         printf("Will print radius normalised by charge\n");
285     }
286
287     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &pbcType, &x, nullptr, box, TRUE);
288     get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
289
290     if (nmol > gnx || gnx % nmol != 0)
291     {
292         gmx_fatal(FARGS, "The number of atoms in the group (%d) is not a multiple of nmol (%d)", gnx, nmol);
293     }
294     nam = gnx / nmol;
295
296     natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
297     snew(x_s, natoms);
298
299     j  = 0;
300     t0 = t;
301     if (bQ)
302     {
303         out = xvgropen(ftp2fn(efXVG, NFILE, fnm),
304                        "Radius of Charge (total and around axes)",
305                        "Time (ps)",
306                        "Rg (nm)",
307                        oenv);
308     }
309     else if (bMOI)
310     {
311         out = xvgropen(ftp2fn(efXVG, NFILE, fnm),
312                        "Moments of inertia (total and around axes)",
313                        "Time (ps)",
314                        "I (a.m.u. nm\\S2\\N)",
315                        oenv);
316     }
317     else
318     {
319         out = xvgropen(ftp2fn(efXVG, NFILE, fnm),
320                        "Radius of gyration (total and around axes)",
321                        "Time (ps)",
322                        "Rg (nm)",
323                        oenv);
324     }
325     if (bMOI)
326     {
327         xvgr_legend(out, NLEG, legI, oenv);
328     }
329     else
330     {
331         if (bRot)
332         {
333             if (output_env_get_print_xvgr_codes(oenv))
334             {
335                 fprintf(out, "@ subtitle \"Axes are principal component axes\"\n");
336             }
337         }
338         xvgr_legend(out, NLEG, leg, oenv);
339     }
340     if (nz == 0)
341     {
342         gpbc = gmx_rmpbc_init(&top.idef, pbcType, natoms);
343     }
344     do
345     {
346         if (nz == 0)
347         {
348             gmx_rmpbc_copy(gpbc, natoms, box, x, x_s);
349         }
350         gyro = 0;
351         clear_rvec(gvec);
352         clear_rvec(gvec1);
353         clear_rvec(d);
354         clear_rvec(d1);
355         for (mol = 0; mol < nmol; mol++)
356         {
357             tm = sub_xcm(nz == 0 ? x_s : x, nam, index + mol * nam, top.atoms.atom, xcm, bQ);
358             if (nz == 0)
359             {
360                 gyro += calc_gyro(
361                         x_s, nam, index + mol * nam, top.atoms.atom, tm, gvec1, d1, bQ, bRot, bMOI, trans);
362             }
363             else
364             {
365                 calc_gyro_z(x, box, nam, index + mol * nam, top.atoms.atom, nz, t, out);
366             }
367             rvec_inc(gvec, gvec1);
368             rvec_inc(d, d1);
369         }
370         if (nmol > 0)
371         {
372             gyro /= nmol;
373             svmul(1.0 / nmol, gvec, gvec);
374             svmul(1.0 / nmol, d, d);
375         }
376
377         if (nz == 0)
378         {
379             if (bRot)
380             {
381                 if (j >= max_moi)
382                 {
383                     max_moi += delta_moi;
384                     for (m = 0; (m < DIM); m++)
385                     {
386                         srenew(moi_trans[m], max_moi * DIM);
387                     }
388                 }
389                 for (m = 0; (m < DIM); m++)
390                 {
391                     copy_rvec(trans[m], moi_trans[m] + DIM * j);
392                 }
393                 fprintf(out, "%10g  %10g  %10g  %10g  %10g\n", t, gyro, d[XX], d[YY], d[ZZ]);
394             }
395             else
396             {
397                 fprintf(out, "%10g  %10g  %10g  %10g  %10g\n", t, gyro, gvec[XX], gvec[YY], gvec[ZZ]);
398             }
399         }
400         j++;
401     } while (read_next_x(oenv, status, &t, x, box));
402     close_trx(status);
403     if (nz == 0)
404     {
405         gmx_rmpbc_done(gpbc);
406     }
407
408     xvgrclose(out);
409
410     if (bACF)
411     {
412         int mode = eacVector;
413
414         do_autocorr(opt2fn("-acf", NFILE, fnm),
415                     oenv,
416                     "Moment of inertia vector ACF",
417                     j,
418                     3,
419                     moi_trans,
420                     (t - t0) / j,
421                     mode,
422                     FALSE);
423         do_view(oenv, opt2fn("-acf", NFILE, fnm), "-nxy");
424     }
425
426     do_view(oenv, ftp2fn(efXVG, NFILE, fnm), "-nxy");
427
428     return 0;
429 }