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