Update instructions in containers.rst
[gromacs.git] / src / gromacs / gmxana / gmx_gyrate.cpp
blob6df19ae6fefeb385c95e1a13455e2c4b519df7d2
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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, 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.
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.
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.
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.
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.
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.
38 #include "gmxpre.h"
40 #include <cmath>
41 #include <cstring>
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/utilities.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/pbcutil/rmpbc.h"
56 #include "gromacs/topology/index.h"
57 #include "gromacs/topology/topology.h"
58 #include "gromacs/utility/arraysize.h"
59 #include "gromacs/utility/fatalerror.h"
60 #include "gromacs/utility/futil.h"
61 #include "gromacs/utility/smalloc.h"
63 static real calc_gyro(rvec x[],
64 int gnx,
65 int index[],
66 t_atom atom[],
67 real tm,
68 rvec gvec,
69 rvec d,
70 gmx_bool bQ,
71 gmx_bool bRot,
72 gmx_bool bMOI,
73 matrix trans)
75 int i, ii, m;
76 real gyro, dx2, m0, Itot;
77 rvec comp;
79 if (bRot)
81 principal_comp(gnx, index, atom, x, trans, d);
82 Itot = norm(d);
83 if (bMOI)
85 return Itot;
87 for (m = 0; (m < DIM); m++)
89 d[m] = std::sqrt(d[m] / tm);
91 /* rotate_atoms(gnx,index,x,trans); */
93 clear_rvec(comp);
94 for (i = 0; (i < gnx); i++)
96 ii = index[i];
97 if (bQ)
99 m0 = std::abs(atom[ii].q);
101 else
103 m0 = atom[ii].m;
105 for (m = 0; (m < DIM); m++)
107 dx2 = x[ii][m] * x[ii][m];
108 comp[m] += dx2 * m0;
111 gyro = comp[XX] + comp[YY] + comp[ZZ];
113 for (m = 0; (m < DIM); m++)
115 gvec[m] = std::sqrt((gyro - comp[m]) / tm);
118 return std::sqrt(gyro / tm);
121 static void calc_gyro_z(rvec x[], matrix box, int gnx, const int index[], t_atom atom[], int nz, real time, FILE* out)
123 static dvec* inertia = nullptr;
124 static double* tm = nullptr;
125 int i, ii, j, zi;
126 real zf, w, sdet, e1, e2;
128 if (inertia == nullptr)
130 snew(inertia, nz);
131 snew(tm, nz);
134 for (i = 0; i < nz; i++)
136 clear_dvec(inertia[i]);
137 tm[i] = 0;
140 for (i = 0; (i < gnx); i++)
142 ii = index[i];
143 zf = nz * x[ii][ZZ] / box[ZZ][ZZ];
144 if (zf >= nz)
146 zf -= nz;
148 if (zf < 0)
150 zf += nz;
152 for (j = 0; j < 2; j++)
154 zi = static_cast<int>(zf + j);
155 if (zi == nz)
157 zi = 0;
159 w = atom[ii].m * (1 + std::cos(M_PI * (zf - zi)));
160 inertia[zi][0] += w * gmx::square(x[ii][YY]);
161 inertia[zi][1] += w * gmx::square(x[ii][XX]);
162 inertia[zi][2] -= w * x[ii][XX] * x[ii][YY];
163 tm[zi] += w;
166 fprintf(out, "%10g", time);
167 for (j = 0; j < nz; j++)
169 for (i = 0; i < 3; i++)
171 inertia[j][i] /= tm[j];
173 sdet = std::sqrt(gmx::square(inertia[j][0] - inertia[j][1]) + 4 * gmx::square(inertia[j][2]));
174 e1 = std::sqrt(0.5 * (inertia[j][0] + inertia[j][1] + sdet));
175 e2 = std::sqrt(0.5 * (inertia[j][0] + inertia[j][1] - sdet));
176 fprintf(out, " %5.3f %5.3f", e1, e2);
178 fprintf(out, "\n");
182 int gmx_gyrate(int argc, char* argv[])
184 const char* desc[] = {
185 "[THISMODULE] computes the radius of gyration of a molecule",
186 "and the radii of gyration about the [IT]x[it]-, [IT]y[it]- and [IT]z[it]-axes,",
187 "as a function of time. The atoms are explicitly mass weighted.[PAR]",
188 "The axis components corresponds to the mass-weighted root-mean-square",
189 "of the radii components orthogonal to each axis, for example:[PAR]",
190 "Rg(x) = sqrt((sum_i m_i (R_i(y)^2 + R_i(z)^2))/(sum_i m_i)).[PAR]",
191 "With the [TT]-nmol[tt] option the radius of gyration will be calculated",
192 "for multiple molecules by splitting the analysis group in equally",
193 "sized parts.[PAR]",
194 "With the option [TT]-nz[tt] 2D radii of gyration in the [IT]x-y[it] plane",
195 "of slices along the [IT]z[it]-axis are calculated."
197 static int nmol = 1, nz = 0;
198 static gmx_bool bQ = FALSE, bRot = FALSE, bMOI = FALSE;
199 t_pargs pa[] = {
200 { "-nmol", FALSE, etINT, { &nmol }, "The number of molecules to analyze" },
201 { "-q",
202 FALSE,
203 etBOOL,
204 { &bQ },
205 "Use absolute value of the charge of an atom as weighting factor instead of mass" },
206 { "-p",
207 FALSE,
208 etBOOL,
209 { &bRot },
210 "Calculate the radii of gyration about the principal axes." },
211 { "-moi",
212 FALSE,
213 etBOOL,
214 { &bMOI },
215 "Calculate the moments of inertia (defined by the principal axes)." },
216 { "-nz",
217 FALSE,
218 etINT,
219 { &nz },
220 "Calculate the 2D radii of gyration of this number of slices along the z-axis" },
222 FILE* out;
223 t_trxstatus* status;
224 t_topology top;
225 PbcType pbcType;
226 rvec * x, *x_s;
227 rvec xcm, gvec, gvec1;
228 matrix box, trans;
229 gmx_bool bACF;
230 real** moi_trans = nullptr;
231 int max_moi = 0, delta_moi = 100;
232 rvec d, d1; /* eigenvalues of inertia tensor */
233 real t, t0, tm, gyro;
234 int natoms;
235 char* grpname;
236 int j, m, gnx, nam, mol;
237 int* index;
238 gmx_output_env_t* oenv;
239 gmx_rmpbc_t gpbc = nullptr;
240 const char* leg[] = { "Rg", "Rg\\sX\\N", "Rg\\sY\\N", "Rg\\sZ\\N" };
241 const char* legI[] = { "Itot", "I1", "I2", "I3" };
242 #define NLEG asize(leg)
243 t_filenm fnm[] = {
244 { efTRX, "-f", nullptr, ffREAD }, { efTPS, nullptr, nullptr, ffREAD },
245 { efNDX, nullptr, nullptr, ffOPTRD }, { efXVG, nullptr, "gyrate", ffWRITE },
246 { efXVG, "-acf", "moi-acf", ffOPTWR },
248 #define NFILE asize(fnm)
249 int npargs;
250 t_pargs* ppa;
252 npargs = asize(pa);
253 ppa = add_acf_pargs(&npargs, pa);
255 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW, NFILE, fnm, npargs, ppa,
256 asize(desc), desc, 0, nullptr, &oenv))
258 sfree(ppa);
259 return 0;
261 bACF = opt2bSet("-acf", NFILE, fnm);
262 if (bACF && nmol != 1)
264 gmx_fatal(FARGS, "Can only do acf with nmol=1");
266 bRot = bRot || bMOI || bACF;
268 if (nz > 0)
269 bMOI = TRUE;
271 if (bRot)
273 printf("Will rotate system along principal axes\n");
274 snew(moi_trans, DIM);
276 if (bMOI)
278 printf("Will print moments of inertia\n");
279 bQ = FALSE;
281 if (bQ)
283 printf("Will print radius normalised by charge\n");
286 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &pbcType, &x, nullptr, box, TRUE);
287 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
289 if (nmol > gnx || gnx % nmol != 0)
291 gmx_fatal(FARGS, "The number of atoms in the group (%d) is not a multiple of nmol (%d)", gnx, nmol);
293 nam = gnx / nmol;
295 natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
296 snew(x_s, natoms);
298 j = 0;
299 t0 = t;
300 if (bQ)
302 out = xvgropen(ftp2fn(efXVG, NFILE, fnm), "Radius of Charge (total and around axes)",
303 "Time (ps)", "Rg (nm)", oenv);
305 else if (bMOI)
307 out = xvgropen(ftp2fn(efXVG, NFILE, fnm), "Moments of inertia (total and around axes)",
308 "Time (ps)", "I (a.m.u. nm\\S2\\N)", oenv);
310 else
312 out = xvgropen(ftp2fn(efXVG, NFILE, fnm), "Radius of gyration (total and around axes)",
313 "Time (ps)", "Rg (nm)", oenv);
315 if (bMOI)
317 xvgr_legend(out, NLEG, legI, oenv);
319 else
321 if (bRot)
323 if (output_env_get_print_xvgr_codes(oenv))
325 fprintf(out, "@ subtitle \"Axes are principal component axes\"\n");
328 xvgr_legend(out, NLEG, leg, oenv);
330 if (nz == 0)
332 gpbc = gmx_rmpbc_init(&top.idef, pbcType, natoms);
336 if (nz == 0)
338 gmx_rmpbc_copy(gpbc, natoms, box, x, x_s);
340 gyro = 0;
341 clear_rvec(gvec);
342 clear_rvec(gvec1);
343 clear_rvec(d);
344 clear_rvec(d1);
345 for (mol = 0; mol < nmol; mol++)
347 tm = sub_xcm(nz == 0 ? x_s : x, nam, index + mol * nam, top.atoms.atom, xcm, bQ);
348 if (nz == 0)
350 gyro += calc_gyro(x_s, nam, index + mol * nam, top.atoms.atom, tm, gvec1, d1, bQ,
351 bRot, bMOI, trans);
353 else
355 calc_gyro_z(x, box, nam, index + mol * nam, top.atoms.atom, nz, t, out);
357 rvec_inc(gvec, gvec1);
358 rvec_inc(d, d1);
360 if (nmol > 0)
362 gyro /= nmol;
363 svmul(1.0 / nmol, gvec, gvec);
364 svmul(1.0 / nmol, d, d);
367 if (nz == 0)
369 if (bRot)
371 if (j >= max_moi)
373 max_moi += delta_moi;
374 for (m = 0; (m < DIM); m++)
376 srenew(moi_trans[m], max_moi * DIM);
379 for (m = 0; (m < DIM); m++)
381 copy_rvec(trans[m], moi_trans[m] + DIM * j);
383 fprintf(out, "%10g %10g %10g %10g %10g\n", t, gyro, d[XX], d[YY], d[ZZ]);
385 else
387 fprintf(out, "%10g %10g %10g %10g %10g\n", t, gyro, gvec[XX], gvec[YY], gvec[ZZ]);
390 j++;
391 } while (read_next_x(oenv, status, &t, x, box));
392 close_trx(status);
393 if (nz == 0)
395 gmx_rmpbc_done(gpbc);
398 xvgrclose(out);
400 if (bACF)
402 int mode = eacVector;
404 do_autocorr(opt2fn("-acf", NFILE, fnm), oenv, "Moment of inertia vector ACF", j, 3,
405 moi_trans, (t - t0) / j, mode, FALSE);
406 do_view(oenv, opt2fn("-acf", NFILE, fnm), "-nxy");
409 do_view(oenv, ftp2fn(efXVG, NFILE, fnm), "-nxy");
411 return 0;