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,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.
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.
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.
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.
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.
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.
44 #include "gromacs/commandline/pargs.h"
45 #include "gromacs/commandline/viewit.h"
46 #include "gromacs/correlationfunctions/autocorr.h"
47 #include "gromacs/fileio/trrio.h"
48 #include "gromacs/fileio/xvgr.h"
49 #include "gromacs/gmxana/gmx_ana.h"
50 #include "gromacs/gmxana/gstat.h"
51 #include "gromacs/math/functions.h"
52 #include "gromacs/math/units.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/topology/index.h"
55 #include "gromacs/utility/arraysize.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/futil.h"
58 #include "gromacs/utility/gmxassert.h"
59 #include "gromacs/utility/pleasecite.h"
60 #include "gromacs/utility/smalloc.h"
62 static void dump_dih_trr(int nframes, int nangles, real** dih, const char* fn, real* time)
64 int i, j, k, l, m, na;
67 matrix box = { { 2, 0, 0 }, { 0, 2, 0 }, { 0, 0, 2 } };
78 printf("There are %d dihedrals. Will fill %d atom positions with cos/sin\n", nangles, na);
80 fio = gmx_trr_open(fn, "w");
81 for (i = 0; (i < nframes); i++)
84 for (j = 0; (j < nangles); j++)
86 for (m = 0; (m < 2); m++)
88 // This is just because the compler and static-analyzer cannot
89 // know that dih[j][i] is always valid. Since it occurs in the innermost
90 // loop over angles and will only trigger on coding errors, we
91 // only enable it for debug builds.
92 GMX_ASSERT(dih != nullptr && dih[j] != nullptr, "Incorrect dihedral array data");
93 x[k][l] = (m == 0) ? std::cos(dih[j][i]) : std::sin(dih[j][i]);
102 gmx_trr_write_frame(fio, i, time[i], 0, box, na, x, nullptr, nullptr);
108 int gmx_g_angle(int argc, char* argv[])
110 static const char* desc[] = {
111 "[THISMODULE] computes the angle distribution for a number of angles",
112 "or dihedrals.[PAR]",
113 "With option [TT]-ov[tt], you can plot the average angle of",
114 "a group of angles as a function of time. With the [TT]-all[tt] option,",
115 "the first graph is the average and the rest are the individual angles.[PAR]",
116 "With the [TT]-of[tt] option, [THISMODULE] also calculates the fraction of trans",
117 "dihedrals (only for dihedrals) as function of time, but this is",
118 "probably only fun for a select few.[PAR]",
119 "With option [TT]-oc[tt], a dihedral correlation function is calculated.[PAR]",
120 "It should be noted that the index file must contain",
121 "atom triplets for angles or atom quadruplets for dihedrals.",
122 "If this is not the case, the program will crash.[PAR]",
123 "With option [TT]-or[tt], a trajectory file is dumped containing cos and",
124 "sin of selected dihedral angles, which subsequently can be used as",
125 "input for a principal components analysis using [gmx-covar].[PAR]",
126 "Option [TT]-ot[tt] plots when transitions occur between",
127 "dihedral rotamers of multiplicity 3 and [TT]-oh[tt]",
128 "records a histogram of the times between such transitions,",
129 "assuming the input trajectory frames are equally spaced in time."
131 static const char* opt[] = { nullptr, "angle", "dihedral", "improper", "ryckaert-bellemans",
133 static gmx_bool bALL = FALSE, bChandler = FALSE, bAverCorr = FALSE, bPBC = TRUE;
134 static real binwidth = 1;
136 { "-type", FALSE, etENUM, { opt }, "Type of angle to analyse" },
141 "Plot all angles separately in the averages file, in the order of appearance in the "
147 "binwidth (degrees) for calculating the distribution" },
148 { "-periodic", FALSE, etBOOL, { &bPBC }, "Print dihedral angles modulo 360 degrees" },
153 "Use Chandler correlation function (N[trans] = 1, N[gauche] = 0) rather than cosine "
154 "correlation function. Trans is defined as phi < -60 or phi > 60." },
159 "Average the correlation functions for the individual angles/dihedrals" }
161 static const char* bugs[] = {
162 "Counting transitions only works for dihedrals with multiplicity 3"
170 real maxang, S2, norm_fac, maxstat;
172 int nframes, maxangstat, mult, *angstat;
173 int i, j, nangles, first, last;
174 gmx_bool bAver, bRb, bPeriodic, bFrac, /* calculate fraction too? */
175 bTrans, /* worry about transtions too? */
176 bCorr; /* correlation function ? */
177 real aver, aver2, aversig; /* fraction trans dihedrals */
180 real** dih = nullptr; /* mega array with all dih. angles at all times*/
181 real * time, *trans_frac, *aver_angle;
182 t_filenm fnm[] = { { efTRX, "-f", nullptr, ffREAD }, { efNDX, nullptr, "angle", ffREAD },
183 { efXVG, "-od", "angdist", ffWRITE }, { efXVG, "-ov", "angaver", ffOPTWR },
184 { efXVG, "-of", "dihfrac", ffOPTWR }, { efXVG, "-ot", "dihtrans", ffOPTWR },
185 { efXVG, "-oh", "trhisto", ffOPTWR }, { efXVG, "-oc", "dihcorr", ffOPTWR },
186 { efTRR, "-or", nullptr, ffOPTWR } };
187 #define NFILE asize(fnm)
190 gmx_output_env_t* oenv;
193 ppa = add_acf_pargs(&npargs, pa);
194 if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME, NFILE, fnm, npargs, ppa,
195 asize(desc), desc, asize(bugs), bugs, &oenv))
205 GMX_RELEASE_ASSERT(opt[0] != nullptr,
206 "Internal option inconsistency; opt[0]==NULL after processing");
216 case 'r': bRb = TRUE; break;
219 if (opt2bSet("-or", NFILE, fnm))
223 gmx_fatal(FARGS, "Can not combine angles with trr dump");
227 please_cite(stdout, "Mu2005a");
231 /* Calculate bin size */
232 maxangstat = gmx::roundToInt(maxang / binwidth);
233 binwidth = maxang / maxangstat;
235 rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
236 nangles = isize / mult;
237 if ((isize % mult) != 0)
240 "number of index elements not multiple of %d, "
241 "these can not be %s\n",
242 mult, (mult == 3) ? "angle triplets" : "dihedral quadruplets");
246 /* Check whether specific analysis has to be performed */
247 bCorr = opt2bSet("-oc", NFILE, fnm);
248 bAver = opt2bSet("-ov", NFILE, fnm);
249 bTrans = opt2bSet("-ot", NFILE, fnm);
250 bFrac = opt2bSet("-of", NFILE, fnm);
251 if (bTrans && opt[0][0] != 'd')
253 fprintf(stderr, "Option -ot should only accompany -type dihedral. Disabling -ot.\n");
257 if (bChandler && !bCorr)
266 " calculating fractions as defined in this program\n"
267 "makes sense for Ryckaert Bellemans dihs. only. Ignoring -of\n\n");
271 if ((bTrans || bFrac || bCorr) && mult == 3)
274 "Can only do transition, fraction or correlation\n"
275 "on dihedrals. Select -d\n");
279 * We need to know the nr of frames so we can allocate memory for an array
280 * with all dihedral angles at all timesteps. Works for me.
282 if (bTrans || bCorr || bALL || opt2bSet("-or", NFILE, fnm))
287 snew(angstat, maxangstat);
289 read_ang_dih(ftp2fn(efTRX, NFILE, fnm), (mult == 3),
290 bALL || bCorr || bTrans || opt2bSet("-or", NFILE, fnm), bRb, bPBC, maxangstat,
291 angstat, &nframes, &time, isize, index, &trans_frac, &aver_angle, dih, oenv);
293 dt = (time[nframes - 1] - time[0]) / (nframes - 1);
297 sprintf(title, "Average Angle: %s", grpname);
298 out = xvgropen(opt2fn("-ov", NFILE, fnm), title, "Time (ps)", "Angle (degrees)", oenv);
299 for (i = 0; (i < nframes); i++)
301 fprintf(out, "%10.5f %8.3f", time[i], aver_angle[i] * RAD2DEG);
304 for (j = 0; (j < nangles); j++)
309 fprintf(out, " %8.3f", std::atan2(std::sin(dd), std::cos(dd)) * RAD2DEG);
313 fprintf(out, " %8.3f", dih[j][i] * RAD2DEG);
321 if (opt2bSet("-or", NFILE, fnm))
323 dump_dih_trr(nframes, nangles, dih, opt2fn("-or", NFILE, fnm), time);
328 sprintf(title, "Trans fraction: %s", grpname);
329 out = xvgropen(opt2fn("-of", NFILE, fnm), title, "Time (ps)", "Fraction", oenv);
331 for (i = 0; (i < nframes); i++)
333 fprintf(out, "%10.5f %10.3f\n", time[i], trans_frac[i]);
334 tfrac += trans_frac[i];
339 fprintf(stderr, "Average trans fraction: %g\n", tfrac);
345 ana_dih_trans(opt2fn("-ot", NFILE, fnm), opt2fn("-oh", NFILE, fnm), dih, nframes, nangles,
346 grpname, time, bRb, oenv);
351 /* Autocorrelation function */
354 fprintf(stderr, "Not enough frames for correlation function\n");
361 real dval, sixty = DEG2RAD * 60;
364 for (i = 0; (i < nangles); i++)
366 for (j = 0; (j < nframes); j++)
371 bTest = (dval > -sixty) && (dval < sixty);
375 bTest = (dval < -sixty) || (dval > sixty);
379 dih[i][j] = dval - tfrac;
396 do_autocorr(opt2fn("-oc", NFILE, fnm), oenv, "Dihedral Autocorrelation Function",
397 nframes, nangles, dih, dt, mode, bAverCorr);
402 /* Determine the non-zero part of the distribution */
403 for (first = 0; (first < maxangstat - 1) && (angstat[first + 1] == 0); first++) {}
404 for (last = maxangstat - 1; (last > 0) && (angstat[last - 1] == 0); last--) {}
407 for (i = 0; (i < nframes); i++)
409 aver += RAD2DEG * aver_angle[i];
410 aver2 += gmx::square(RAD2DEG * aver_angle[i]);
414 aversig = std::sqrt(aver2 - gmx::square(aver));
415 printf("Found points in the range from %d to %d (max %d)\n", first, last, maxangstat);
416 printf(" < angle > = %g\n", aver);
417 printf("< angle^2 > = %g\n", aver2);
418 printf("Std. Dev. = %g\n", aversig);
422 sprintf(title, "Angle Distribution: %s", grpname);
426 sprintf(title, "Dihedral Distribution: %s", grpname);
428 calc_distribution_props(maxangstat, angstat, -180.0, 0, nullptr, &S2);
429 fprintf(stderr, "Order parameter S^2 = %g\n", S2);
432 bPeriodic = (mult == 4) && (first == 0) && (last == maxangstat - 1);
434 out = xvgropen(opt2fn("-od", NFILE, fnm), title, "Degrees", "", oenv);
435 if (output_env_get_print_xvgr_codes(oenv))
437 fprintf(out, "@ subtitle \"average angle: %g\\So\\N\"\n", aver);
439 norm_fac = 1.0 / (nangles * nframes * binwidth);
443 for (i = first; (i <= last); i++)
445 maxstat = std::max(maxstat, angstat[i] * norm_fac);
447 if (output_env_get_print_xvgr_codes(oenv))
449 fprintf(out, "@with g0\n");
450 fprintf(out, "@ world xmin -180\n");
451 fprintf(out, "@ world xmax 180\n");
452 fprintf(out, "@ world ymin 0\n");
453 fprintf(out, "@ world ymax %g\n", maxstat * 1.05);
454 fprintf(out, "@ xaxis tick major 60\n");
455 fprintf(out, "@ xaxis tick minor 30\n");
456 fprintf(out, "@ yaxis tick major 0.005\n");
457 fprintf(out, "@ yaxis tick minor 0.0025\n");
460 for (i = first; (i <= last); i++)
462 fprintf(out, "%10g %10f\n", i * binwidth + 180.0 - maxang, angstat[i] * norm_fac);
466 /* print first bin again as last one */
467 fprintf(out, "%10g %10f\n", 180.0, angstat[0] * norm_fac);
472 do_view(oenv, opt2fn("-od", NFILE, fnm), "-nxy");
475 do_view(oenv, opt2fn("-ov", NFILE, fnm), "-nxy");