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.
45 #include "gromacs/fileio/confio.h"
46 #include "gromacs/fileio/trxio.h"
47 #include "gromacs/fileio/xvgr.h"
48 #include "gromacs/gmxana/angle_correction.h"
49 #include "gromacs/gmxana/gstat.h"
50 #include "gromacs/listed_forces/bonded.h"
51 #include "gromacs/math/functions.h"
52 #include "gromacs/math/units.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/math/vecdump.h"
55 #include "gromacs/pbcutil/pbc.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/smalloc.h"
59 void print_one(const gmx_output_env_t* oenv,
69 char buf[256], t2[256];
72 sprintf(buf, "%s%s.xvg", base, name);
73 fprintf(stderr, "\rPrinting %s ", buf);
75 sprintf(t2, "%s %s", title, name);
76 fp = xvgropen(buf, t2, "Time (ps)", ylabel, oenv);
77 for (k = 0; (k < nf); k++)
79 fprintf(fp, "%10g %10g\n", time[k], data[k]);
84 static int calc_RBbin(real phi, int gmx_unused multiplicity, real gmx_unused core_frac)
86 /* multiplicity and core_frac NOT used,
87 * just given to enable use of pt-to-fn in caller low_ana_dih_trans*/
88 static const real r30 = M_PI / 6.0;
89 static const real r90 = M_PI / 2.0;
90 static const real r150 = M_PI * 5.0 / 6.0;
92 if ((phi < r30) && (phi > -r30))
96 else if ((phi > -r150) && (phi < -r90))
100 else if ((phi < r150) && (phi > r90))
107 static int calc_Nbin(real phi, int multiplicity, real core_frac)
109 static const real r360 = 360 * DEG2RAD;
110 real rot_width, core_width, core_offset, low, hi;
112 /* with multiplicity 3 and core_frac 0.5
113 * 0<g(-)<120, 120<t<240, 240<g(+)<360
114 * 0< bin0 < 30, 30<bin1<90, 90<bin0<150, 150<bin2<210, 210<bin0<270, 270<bin3<330, 330<bin0<360
115 * so with multiplicity 3, bin1 is core g(-), bin2 is core t, bin3 is
116 core g(+), bin0 is between rotamers */
122 rot_width = 360. / multiplicity;
123 core_width = core_frac * rot_width;
124 core_offset = (rot_width - core_width) / 2.0;
125 for (bin = 1; bin <= multiplicity; bin++)
127 low = ((bin - 1) * rot_width) + core_offset;
128 hi = ((bin - 1) * rot_width) + core_offset + core_width;
131 if ((phi > low) && (phi < hi))
139 void ana_dih_trans(const char* fn_trans,
140 const char* fn_histo,
147 const gmx_output_env_t* oenv)
149 /* just a wrapper; declare extra args, then chuck away at end. */
157 snew(multiplicity, nangles);
158 for (k = 0; (k < nangles); k++)
163 low_ana_dih_trans(TRUE, fn_trans, TRUE, fn_histo, maxchi, dih, nlist, dlist, nframes, nangles,
164 grpname, multiplicity, time, bRb, 0.5, oenv);
169 void low_ana_dih_trans(gmx_bool bTrans,
170 const char* fn_trans,
172 const char* fn_histo,
184 const gmx_output_env_t* oenv)
189 int i, j, k, Dih, ntrans;
190 int cur_bin, new_bin;
193 int (*calc_bin)(real, int, real);
200 /* Assumes the frames are equally spaced in time */
201 dt = (time[nframes - 1] - time[0]) / (nframes - 1);
203 /* Analysis of dihedral transitions */
204 fprintf(stderr, "Now calculating transitions...\n");
208 calc_bin = calc_RBbin;
212 calc_bin = calc_Nbin;
215 for (k = 0; k < NROT; k++)
217 snew(rot_occ[k], nangles);
218 for (i = 0; (i < nangles); i++)
226 /* dih[i][j] is the dihedral angle i in frame j */
228 for (i = 0; (i < nangles); i++)
233 mind = maxd = prev = dih[i][0];
235 cur_bin = calc_bin(dih[i][0], multiplicity[i], core_frac);
236 rot_occ[cur_bin][i]++;
238 for (j = 1; (j < nframes); j++)
240 new_bin = calc_bin(dih[i][j], multiplicity[i], core_frac);
241 rot_occ[new_bin][i]++;
247 else if ((new_bin != 0) && (cur_bin != new_bin))
255 /* why is all this md rubbish periodic? Remove 360 degree periodicity */
256 if ((dih[i][j] - prev) > M_PI)
258 dih[i][j] -= 2 * M_PI;
260 else if ((dih[i][j] - prev) < -M_PI)
262 dih[i][j] += 2 * M_PI;
267 mind = std::min(mind, dih[i][j]);
268 maxd = std::max(maxd, dih[i][j]);
269 if ((maxd - mind) > 2 * M_PI / 3) /* or 120 degrees, assuming */
270 { /* multiplicity 3. Not so general.*/
273 maxd = mind = dih[i][j]; /* get ready for next transition */
278 for (k = 0; k < NROT; k++)
280 rot_occ[k][i] /= nframes;
283 fprintf(stderr, "Total number of transitions: %10d\n", ntrans);
286 ttime = (dt * nframes * nangles) / ntrans;
287 fprintf(stderr, "Time between transitions: %10.3f ps\n", ttime);
290 /* new by grs - copy transitions from tr_h[] to dlist->ntr[]
291 * and rotamer populations from rot_occ to dlist->rot_occ[]
292 * based on fn histogramming in g_chi. diff roles for i and j here */
295 for (Dih = 0; (Dih < NONCHI + maxchi); Dih++)
297 for (i = 0; (i < nlist); i++)
299 if (((Dih < edOmega)) || ((Dih == edOmega) && (has_dihedral(edOmega, &(dlist[i]))))
300 || ((Dih > edOmega) && (dlist[i].atm.Cn[Dih - NONCHI + 3] != -1)))
302 /* grs debug printf("Not OK? i %d j %d Dih %d \n", i, j, Dih) ; */
303 dlist[i].ntr[Dih] = tr_h[j];
304 for (k = 0; k < NROT; k++)
306 dlist[i].rot_occ[Dih][k] = rot_occ[k][j];
313 /* end addition by grs */
317 sprintf(title, "Number of transitions: %s", grpname);
318 fp = xvgropen(fn_trans, title, "Time (ps)", "# transitions/timeframe", oenv);
319 for (j = 0; (j < nframes); j++)
321 fprintf(fp, "%10.3f %10d\n", time[j], tr_f[j]);
326 /* Compute histogram from # transitions per dihedral */
328 for (j = 0; (j < nframes); j++)
332 for (i = 0; (i < nangles); i++)
336 for (j = nframes; ((tr_f[j - 1] == 0) && (j > 0)); j--) {}
338 ttime = dt * nframes;
341 sprintf(title, "Transition time: %s", grpname);
342 fp = xvgropen(fn_histo, title, "Time (ps)", "#", oenv);
343 for (i = j - 1; (i > 0); i--)
347 fprintf(fp, "%10.3f %10d\n", ttime / i, tr_f[i]);
355 for (k = 0; k < NROT; k++)
361 void mk_multiplicity_lookup(int* multiplicity, int maxchi, int nlist, t_dlist dlist[], int nangles)
363 /* new by grs - for dihedral j (as in dih[j]) get multiplicity from dlist
364 * and store in multiplicity[j]
371 for (Dih = 0; (Dih < NONCHI + maxchi); Dih++)
373 for (i = 0; (i < nlist); i++)
375 std::strncpy(name, dlist[i].name, 3);
377 if (((Dih < edOmega)) || ((Dih == edOmega) && (has_dihedral(edOmega, &(dlist[i]))))
378 || ((Dih > edOmega) && (dlist[i].atm.Cn[Dih - NONCHI + 3] != -1)))
380 /* default - we will correct the rest below */
383 /* make omegas 2fold, though doesn't make much more sense than 3 */
384 if (Dih == edOmega && (has_dihedral(edOmega, &(dlist[i]))))
389 /* dihedrals to aromatic rings, COO, CONH2 or guanidinium are 2fold*/
390 if (Dih > edOmega && (dlist[i].atm.Cn[Dih - NONCHI + 3] != -1))
392 if (((std::strstr(name, "PHE") != nullptr) && (Dih == edChi2))
393 || ((std::strstr(name, "TYR") != nullptr) && (Dih == edChi2))
394 || ((std::strstr(name, "PTR") != nullptr) && (Dih == edChi2))
395 || ((std::strstr(name, "TRP") != nullptr) && (Dih == edChi2))
396 || ((std::strstr(name, "HIS") != nullptr) && (Dih == edChi2))
397 || ((std::strstr(name, "GLU") != nullptr) && (Dih == edChi3))
398 || ((std::strstr(name, "ASP") != nullptr) && (Dih == edChi2))
399 || ((std::strstr(name, "GLN") != nullptr) && (Dih == edChi3))
400 || ((std::strstr(name, "ASN") != nullptr) && (Dih == edChi2))
401 || ((std::strstr(name, "ARG") != nullptr) && (Dih == edChi4)))
412 fprintf(stderr, "WARNING: not all dihedrals found in topology (only %d out of %d)!\n", j, nangles);
414 /* Check for remaining dihedrals */
415 for (; (j < nangles); j++)
421 void mk_chi_lookup(int** lookup, int maxchi, int nlist, t_dlist dlist[])
424 /* by grs. should rewrite everything to use this. (but haven't,
425 * and at mmt only used in get_chi_product_traj
426 * returns the dihed number given the residue number (from-0)
427 * and chi (from-0) nr. -1 for chi undefined for that res (eg gly, ala..)*/
432 /* NONCHI points to chi1, therefore we have to start counting there. */
433 for (Dih = NONCHI; (Dih < NONCHI + maxchi); Dih++)
435 for (i = 0; (i < nlist); i++)
438 if (((Dih < edOmega)) || ((Dih == edOmega) && (has_dihedral(edOmega, &(dlist[i]))))
439 || ((Dih > edOmega) && (dlist[i].atm.Cn[Dih - NONCHI + 3] != -1)))
441 /* grs debug printf("Not OK? i %d j %d Dih %d \n", i, j, Dih) ; */
457 void get_chi_product_traj(real** dih,
470 const gmx_output_env_t* oenv)
473 gmx_bool bRotZero, bHaveChi = FALSE;
474 int accum = 0, index, i, j, k, Xi, n, b;
479 char hisfile[256], histitle[256], *namept;
481 int (*calc_bin)(real, int, real);
483 /* Analysis of dihedral transitions */
484 fprintf(stderr, "Now calculating Chi product trajectories...\n");
488 calc_bin = calc_RBbin;
492 calc_bin = calc_Nbin;
495 snew(chi_prtrj, nframes);
497 /* file for info on all residues */
500 fpall = xvgropen(fnall, "Cumulative Rotamers", "Residue", "Probability", oenv);
504 fpall = xvgropen(fnall, "Cumulative Rotamers", "Residue", "# Counts", oenv);
507 for (i = 0; (i < nlist); i++)
510 /* get nbin, the nr. of cumulative rotamers that need to be considered */
512 for (Xi = 0; Xi < maxchi; Xi++)
514 index = lookup[i][Xi]; /* chi_(Xi+1) of res i (-1 if off end) */
517 n = multiplicity[index];
521 nbin += 1; /* for the "zero rotamer", outside the core region */
523 for (j = 0; (j < nframes); j++)
528 index = lookup[i][0]; /* index into dih of chi1 of res i */
536 b = calc_bin(dih[index][j], multiplicity[index], core_frac);
542 for (Xi = 1; Xi < maxchi; Xi++)
544 index = lookup[i][Xi]; /* chi_(Xi+1) of res i (-1 if off end) */
547 n = multiplicity[index];
548 b = calc_bin(dih[index][j], n, core_frac);
549 accum = n * accum + b - 1;
564 chi_prtrj[j] = accum;
565 if (accum + 1 > nbin)
576 /* print cuml rotamer vs time */
577 print_one(oenv, "chiproduct", dlist[i].name, "chi product for",
578 "cumulative rotamer", nframes, time, chi_prtrj);
581 /* make a histogram pf culm. rotamer occupancy too */
582 snew(chi_prhist, nbin);
583 make_histo(nullptr, nframes, chi_prtrj, nbin, chi_prhist, 0, nbin);
586 sprintf(hisfile, "histo-chiprod%s.xvg", dlist[i].name);
587 sprintf(histitle, "cumulative rotamer distribution for %s", dlist[i].name);
588 fprintf(stderr, " and %s ", hisfile);
589 fp = xvgropen(hisfile, histitle, "number", "", oenv);
590 if (output_env_get_print_xvgr_codes(oenv))
592 fprintf(fp, "@ xaxis tick on\n");
593 fprintf(fp, "@ xaxis tick major 1\n");
594 fprintf(fp, "@ type xy\n");
596 for (k = 0; (k < nbin); k++)
600 fprintf(fp, "%5d %10g\n", k, (1.0 * chi_prhist[k]) / nframes);
604 fprintf(fp, "%5d %10d\n", k, chi_prhist[k]);
607 fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
611 /* and finally print out occupancies to a single file */
612 /* get the gmx from-1 res nr by setting a ptr to the number part
613 * of dlist[i].name - potential bug for 4-letter res names... */
614 namept = dlist[i].name + 3;
615 fprintf(fpall, "%5s ", namept);
616 for (k = 0; (k < nbin); k++)
620 fprintf(fpall, " %10g", (1.0 * chi_prhist[k]) / nframes);
624 fprintf(fpall, " %10d", chi_prhist[k]);
627 fprintf(fpall, "\n");
636 fprintf(stderr, "\n");
639 void calc_distribution_props(int nh, const int histo[], real start, int nkkk, t_karplus kkk[], real* S2)
641 real d, dc, ds, c1, c2, tdc, tds;
642 real fac, ang, invth, Jc;
647 gmx_fatal(FARGS, "No points in histogram (%s, %d)", __FILE__, __LINE__);
651 /* Compute normalisation factor */
653 for (j = 0; (j < nh); j++)
659 for (i = 0; (i < nkkk); i++)
666 for (j = 0; (j < nh); j++)
668 d = invth * histo[j];
669 ang = j * fac - start;
672 ds = d * std::sin(ang);
675 for (i = 0; (i < nkkk); i++)
677 c1 = std::cos(ang + kkk[i].offset);
679 Jc = (kkk[i].A * c2 + kkk[i].B * c1 + kkk[i].C);
680 kkk[i].Jc += histo[j] * Jc;
681 kkk[i].Jcsig += histo[j] * gmx::square(Jc);
684 for (i = 0; (i < nkkk); i++)
687 kkk[i].Jcsig = std::sqrt(kkk[i].Jcsig / th - gmx::square(kkk[i].Jc));
689 *S2 = tdc * tdc + tds * tds;
692 static void calc_angles(struct t_pbc* pbc, int n3, int index[], real ang[], rvec x_s[])
698 for (i = ix = 0; (ix < n3); i++, ix += 3)
700 ang[i] = bond_angle(x_s[index[ix]], x_s[index[ix + 1]], x_s[index[ix + 2]], pbc, r_ij, r_kj,
705 fprintf(debug, "Angle[0]=%g, costh=%g, index0 = %d, %d, %d\n", ang[0], costh, index[0],
707 pr_rvec(debug, 0, "rij", r_ij, DIM, TRUE);
708 pr_rvec(debug, 0, "rkj", r_kj, DIM, TRUE);
712 static real calc_fraction(const real angles[], int nangles)
715 real trans = 0, gauche = 0;
718 for (i = 0; i < nangles; i++)
720 angle = angles[i] * RAD2DEG;
722 if (angle > 135 && angle < 225)
726 else if (angle > 270 && angle < 330)
730 else if (angle < 90 && angle > 30)
735 if (trans + gauche > 0)
737 return trans / (trans + gauche);
745 static void calc_dihs(struct t_pbc* pbc, int n4, const int index[], real ang[], rvec x_s[])
747 int i, ix, t1, t2, t3;
748 rvec r_ij, r_kj, r_kl, m, n;
751 for (i = ix = 0; (ix < n4); i++, ix += 4)
753 aaa = dih_angle(x_s[index[ix]], x_s[index[ix + 1]], x_s[index[ix + 2]], x_s[index[ix + 3]],
754 pbc, r_ij, r_kj, r_kl, m, n, &t1, &t2, &t3);
756 ang[i] = aaa; /* not taking into account ryckaert bellemans yet */
760 void make_histo(FILE* log, int ndata, real data[], int npoints, int histo[], real minx, real maxx)
767 minx = maxx = data[0];
768 for (i = 1; (i < ndata); i++)
770 minx = std::min(minx, data[i]);
771 maxx = std::max(maxx, data[i]);
773 fprintf(log, "Min data: %10g Max data: %10g\n", minx, maxx);
775 dx = npoints / (maxx - minx);
778 fprintf(debug, "Histogramming: ndata=%d, nhisto=%d, minx=%g,maxx=%g,dx=%g\n", ndata,
779 npoints, minx, maxx, dx);
781 for (i = 0; (i < ndata); i++)
783 ind = static_cast<int>((data[i] - minx) * dx);
784 if ((ind >= 0) && (ind < npoints))
790 fprintf(log, "index = %d, data[%d] = %g\n", ind, i, data[i]);
795 void normalize_histo(int npoints, const int histo[], real dx, real normhisto[])
801 for (i = 0; (i < npoints); i++)
807 fprintf(stderr, "Empty histogram!\n");
811 for (i = 0; (i < npoints); i++)
813 normhisto[i] = fac * histo[i];
817 void read_ang_dih(const char* trj_fn,
831 const gmx_output_env_t* oenv)
835 int i, angind, total, teller;
836 int nangles, n_alloc;
837 real t, fraction, pifac, angle;
842 #define prev (1 - cur)
845 read_first_x(oenv, &status, trj_fn, &t, &x, box);
857 snew(angles[cur], nangles);
858 snew(angles[prev], nangles);
860 /* Start the loop over frames */
865 *trans_frac = nullptr;
866 *aver_angle = nullptr;
870 if (teller >= n_alloc)
875 for (i = 0; (i < nangles); i++)
877 srenew(dih[i], n_alloc);
880 srenew(*time, n_alloc);
881 srenew(*trans_frac, n_alloc);
882 srenew(*aver_angle, n_alloc);
889 set_pbc(pbc, -1, box);
894 calc_angles(pbc, isize, index, angles[cur], x);
898 calc_dihs(pbc, isize, index, angles[cur], x);
901 fraction = calc_fraction(angles[cur], nangles);
902 (*trans_frac)[teller] = fraction;
904 /* Change Ryckaert-Bellemans dihedrals to polymer convention
905 * Modified 990913 by Erik:
906 * We actually shouldn't change the convention, since it's
907 * calculated from polymer above, but we change the intervall
908 * from [-180,180] to [0,360].
912 for (i = 0; (i < nangles); i++)
914 if (angles[cur][i] <= 0.0)
916 angles[cur][i] += 2 * M_PI;
921 /* Periodicity in dihedral space... */
924 for (i = 0; (i < nangles); i++)
926 real dd = angles[cur][i];
927 angles[cur][i] = std::atan2(std::sin(dd), std::cos(dd));
934 for (i = 0; (i < nangles); i++)
936 while (angles[cur][i] <= angles[prev][i] - M_PI)
938 angles[cur][i] += 2 * M_PI;
940 while (angles[cur][i] > angles[prev][i] + M_PI)
942 angles[cur][i] -= 2 * M_PI;
951 for (i = 0; (i < nangles); i++)
953 if (!bAngles && i > 0)
955 real diffa = angles[cur][i] - angles[cur][i - 1];
956 diffa = correctRadianAngleRange(diffa);
957 angles[cur][i] = angles[cur][i - 1] + diffa;
960 aa = aa + angles[cur][i];
962 /* angle in rad / 2Pi * max determines bin. bins go from 0 to maxangstat,
963 even though scale goes from -pi to pi (dihedral) or -pi/2 to pi/2
964 (angle) Basically: translate the x-axis by Pi. Translate it back by
968 angle = angles[cur][i];
971 angle = correctRadianAngleRange(angle);
975 /* Update the distribution histogram */
976 angind = gmx::roundToInt((angle * maxangstat) / pifac);
977 if (angind == maxangstat)
981 if ((angind < 0) || (angind >= maxangstat))
983 /* this will never happen */
984 gmx_fatal(FARGS, "angle (%f) index out of range (0..%d) : %d\n", angle, maxangstat, angind);
988 if (angind == maxangstat)
990 fprintf(stderr, "angle %d fr %d = %g\n", i, cur, angle);
996 /* average over all angles */
997 aa = correctRadianAngleRange(aa / nangles);
998 (*aver_angle)[teller] = (aa);
1000 /* this copies all current dih. angles to dih[i], teller is frame */
1003 for (i = 0; i < nangles; i++)
1007 dih[i][teller] = correctRadianAngleRange(angles[cur][i]);
1011 dih[i][teller] = angles[cur][i];
1019 /* Increment loop counter */
1021 } while (read_next_x(oenv, status, &t, x, box));
1026 sfree(angles[prev]);