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, 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.
43 #include "gromacs/bonded/bonded.h"
44 #include "gromacs/fileio/confio.h"
45 #include "gromacs/fileio/trxio.h"
46 #include "gromacs/fileio/xvgr.h"
47 #include "gromacs/gmxana/gstat.h"
48 #include "gromacs/legacyheaders/macros.h"
49 #include "gromacs/legacyheaders/txtdump.h"
50 #include "gromacs/legacyheaders/typedefs.h"
51 #include "gromacs/math/units.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/pbcutil/pbc.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/smalloc.h"
57 void print_one(const output_env_t oenv, const char *base, const char *name,
58 const char *title, const char *ylabel, int nf, real time[],
62 char buf[256], t2[256];
65 sprintf(buf, "%s%s.xvg", base, name);
66 fprintf(stderr, "\rPrinting %s ", buf);
67 sprintf(t2, "%s %s", title, name);
68 fp = xvgropen(buf, t2, "Time (ps)", ylabel, oenv);
69 for (k = 0; (k < nf); k++)
71 fprintf(fp, "%10g %10g\n", time[k], data[k]);
76 static int calc_RBbin(real phi, int gmx_unused multiplicity, real gmx_unused core_frac)
78 /* multiplicity and core_frac NOT used,
79 * just given to enable use of pt-to-fn in caller low_ana_dih_trans*/
80 static const real r30 = M_PI/6.0;
81 static const real r90 = M_PI/2.0;
82 static const real r150 = M_PI*5.0/6.0;
84 if ((phi < r30) && (phi > -r30))
88 else if ((phi > -r150) && (phi < -r90))
92 else if ((phi < r150) && (phi > r90))
99 static int calc_Nbin(real phi, int multiplicity, real core_frac)
101 static const real r360 = 360*DEG2RAD;
102 real rot_width, core_width, core_offset, low, hi;
104 /* with multiplicity 3 and core_frac 0.5
105 * 0<g(-)<120, 120<t<240, 240<g(+)<360
106 * 0< bin0 < 30, 30<bin1<90, 90<bin0<150, 150<bin2<210, 210<bin0<270, 270<bin3<330, 330<bin0<360
107 * so with multiplicity 3, bin1 is core g(-), bin2 is core t, bin3 is
108 core g(+), bin0 is between rotamers */
114 rot_width = 360/multiplicity;
115 core_width = core_frac * rot_width;
116 core_offset = (rot_width - core_width)/2.0;
117 for (bin = 1; bin <= multiplicity; bin++)
119 low = ((bin - 1) * rot_width ) + core_offset;
120 hi = ((bin - 1) * rot_width ) + core_offset + core_width;
123 if ((phi > low) && (phi < hi))
131 void ana_dih_trans(const char *fn_trans, const char *fn_histo,
132 real **dih, int nframes, int nangles,
133 const char *grpname, real *time, gmx_bool bRb,
134 const output_env_t oenv)
136 /* just a wrapper; declare extra args, then chuck away at end. */
144 snew(multiplicity, nangles);
145 for (k = 0; (k < nangles); k++)
150 low_ana_dih_trans(TRUE, fn_trans, TRUE, fn_histo, maxchi,
151 dih, nlist, dlist, nframes,
152 nangles, grpname, multiplicity, time, bRb, 0.5, oenv);
158 void low_ana_dih_trans(gmx_bool bTrans, const char *fn_trans,
159 gmx_bool bHisto, const char *fn_histo, int maxchi,
160 real **dih, int nlist, t_dlist dlist[], int nframes,
161 int nangles, const char *grpname, int multiplicity[],
162 real *time, gmx_bool bRb, real core_frac,
163 const output_env_t oenv)
168 int i, j, k, Dih, ntrans;
169 int cur_bin, new_bin;
172 int (*calc_bin)(real, int, real);
179 /* Assumes the frames are equally spaced in time */
180 dt = (time[nframes-1]-time[0])/(nframes-1);
182 /* Analysis of dihedral transitions */
183 fprintf(stderr, "Now calculating transitions...\n");
187 calc_bin = calc_RBbin;
191 calc_bin = calc_Nbin;
194 for (k = 0; k < NROT; k++)
196 snew(rot_occ[k], nangles);
197 for (i = 0; (i < nangles); i++)
205 /* dih[i][j] is the dihedral angle i in frame j */
207 for (i = 0; (i < nangles); i++)
212 mind = maxd = prev = dih[i][0];
214 cur_bin = calc_bin(dih[i][0], multiplicity[i], core_frac);
215 rot_occ[cur_bin][i]++;
217 for (j = 1; (j < nframes); j++)
219 new_bin = calc_bin(dih[i][j], multiplicity[i], core_frac);
220 rot_occ[new_bin][i]++;
226 else if ((new_bin != 0) && (cur_bin != new_bin))
234 /* why is all this md rubbish periodic? Remove 360 degree periodicity */
235 if ( (dih[i][j] - prev) > M_PI)
239 else if ( (dih[i][j] - prev) < -M_PI)
246 mind = min(mind, dih[i][j]);
247 maxd = max(maxd, dih[i][j]);
248 if ( (maxd - mind) > 2*M_PI/3) /* or 120 degrees, assuming */
249 { /* multiplicity 3. Not so general.*/
252 maxd = mind = dih[i][j]; /* get ready for next transition */
257 for (k = 0; k < NROT; k++)
259 rot_occ[k][i] /= nframes;
262 fprintf(stderr, "Total number of transitions: %10d\n", ntrans);
265 ttime = (dt*nframes*nangles)/ntrans;
266 fprintf(stderr, "Time between transitions: %10.3f ps\n", ttime);
269 /* new by grs - copy transitions from tr_h[] to dlist->ntr[]
270 * and rotamer populations from rot_occ to dlist->rot_occ[]
271 * based on fn histogramming in g_chi. diff roles for i and j here */
274 for (Dih = 0; (Dih < NONCHI+maxchi); Dih++)
276 for (i = 0; (i < nlist); i++)
278 if (((Dih < edOmega) ) ||
279 ((Dih == edOmega) && (has_dihedral(edOmega, &(dlist[i])))) ||
280 ((Dih > edOmega) && (dlist[i].atm.Cn[Dih-NONCHI+3] != -1)))
282 /* grs debug printf("Not OK? i %d j %d Dih %d \n", i, j, Dih) ; */
283 dlist[i].ntr[Dih] = tr_h[j];
284 for (k = 0; k < NROT; k++)
286 dlist[i].rot_occ[Dih][k] = rot_occ[k][j];
293 /* end addition by grs */
297 sprintf(title, "Number of transitions: %s", grpname);
298 fp = xvgropen(fn_trans, title, "Time (ps)", "# transitions/timeframe", oenv);
299 for (j = 0; (j < nframes); j++)
301 fprintf(fp, "%10.3f %10d\n", time[j], tr_f[j]);
306 /* Compute histogram from # transitions per dihedral */
308 for (j = 0; (j < nframes); j++)
312 for (i = 0; (i < nangles); i++)
316 for (j = nframes; ((tr_f[j-1] == 0) && (j > 0)); j--)
324 sprintf(title, "Transition time: %s", grpname);
325 fp = xvgropen(fn_histo, title, "Time (ps)", "#", oenv);
326 for (i = j-1; (i > 0); i--)
330 fprintf(fp, "%10.3f %10d\n", ttime/i, tr_f[i]);
338 for (k = 0; k < NROT; k++)
345 void mk_multiplicity_lookup (int *multiplicity, int maxchi,
346 int nlist, t_dlist dlist[], int nangles)
348 /* new by grs - for dihedral j (as in dih[j]) get multiplicity from dlist
349 * and store in multiplicity[j]
356 for (Dih = 0; (Dih < NONCHI+maxchi); Dih++)
358 for (i = 0; (i < nlist); i++)
360 strncpy(name, dlist[i].name, 3);
362 if (((Dih < edOmega) ) ||
363 ((Dih == edOmega) && (has_dihedral(edOmega, &(dlist[i])))) ||
364 ((Dih > edOmega) && (dlist[i].atm.Cn[Dih-NONCHI+3] != -1)))
366 /* default - we will correct the rest below */
369 /* make omegas 2fold, though doesn't make much more sense than 3 */
370 if (Dih == edOmega && (has_dihedral(edOmega, &(dlist[i]))))
375 /* dihedrals to aromatic rings, COO, CONH2 or guanidinium are 2fold*/
376 if (Dih > edOmega && (dlist[i].atm.Cn[Dih-NONCHI+3] != -1))
378 if ( ((strstr(name, "PHE") != NULL) && (Dih == edChi2)) ||
379 ((strstr(name, "TYR") != NULL) && (Dih == edChi2)) ||
380 ((strstr(name, "PTR") != NULL) && (Dih == edChi2)) ||
381 ((strstr(name, "TRP") != NULL) && (Dih == edChi2)) ||
382 ((strstr(name, "HIS") != NULL) && (Dih == edChi2)) ||
383 ((strstr(name, "GLU") != NULL) && (Dih == edChi3)) ||
384 ((strstr(name, "ASP") != NULL) && (Dih == edChi2)) ||
385 ((strstr(name, "GLN") != NULL) && (Dih == edChi3)) ||
386 ((strstr(name, "ASN") != NULL) && (Dih == edChi2)) ||
387 ((strstr(name, "ARG") != NULL) && (Dih == edChi4)) )
398 fprintf(stderr, "WARNING: not all dihedrals found in topology (only %d out of %d)!\n",
401 /* Check for remaining dihedrals */
402 for (; (j < nangles); j++)
409 void mk_chi_lookup (int **lookup, int maxchi,
410 int nlist, t_dlist dlist[])
413 /* by grs. should rewrite everything to use this. (but haven't,
414 * and at mmt only used in get_chi_product_traj
415 * returns the dihed number given the residue number (from-0)
416 * and chi (from-0) nr. -1 for chi undefined for that res (eg gly, ala..)*/
421 /* NONCHI points to chi1, therefore we have to start counting there. */
422 for (Dih = NONCHI; (Dih < NONCHI+maxchi); Dih++)
424 for (i = 0; (i < nlist); i++)
427 if (((Dih < edOmega) ) ||
428 ((Dih == edOmega) && (has_dihedral(edOmega, &(dlist[i])))) ||
429 ((Dih > edOmega) && (dlist[i].atm.Cn[Dih-NONCHI+3] != -1)))
431 /* grs debug printf("Not OK? i %d j %d Dih %d \n", i, j, Dih) ; */
448 void get_chi_product_traj (real **dih, int nframes, int nlist,
449 int maxchi, t_dlist dlist[], real time[],
450 int **lookup, int *multiplicity, gmx_bool bRb, gmx_bool bNormalize,
451 real core_frac, gmx_bool bAll, const char *fnall,
452 const output_env_t oenv)
455 gmx_bool bRotZero, bHaveChi = FALSE;
456 int accum = 0, index, i, j, k, Xi, n, b;
461 char hisfile[256], histitle[256], *namept;
463 int (*calc_bin)(real, int, real);
465 /* Analysis of dihedral transitions */
466 fprintf(stderr, "Now calculating Chi product trajectories...\n");
470 calc_bin = calc_RBbin;
474 calc_bin = calc_Nbin;
477 snew(chi_prtrj, nframes);
479 /* file for info on all residues */
482 fpall = xvgropen(fnall, "Cumulative Rotamers", "Residue", "Probability", oenv);
486 fpall = xvgropen(fnall, "Cumulative Rotamers", "Residue", "# Counts", oenv);
489 for (i = 0; (i < nlist); i++)
492 /* get nbin, the nr. of cumulative rotamers that need to be considered */
494 for (Xi = 0; Xi < maxchi; Xi++)
496 index = lookup[i][Xi]; /* chi_(Xi+1) of res i (-1 if off end) */
499 n = multiplicity[index];
503 nbin += 1; /* for the "zero rotamer", outside the core region */
505 for (j = 0; (j < nframes); j++)
510 index = lookup[i][0]; /* index into dih of chi1 of res i */
519 b = calc_bin(dih[index][j], multiplicity[index], core_frac);
525 for (Xi = 1; Xi < maxchi; Xi++)
527 index = lookup[i][Xi]; /* chi_(Xi+1) of res i (-1 if off end) */
530 n = multiplicity[index];
531 b = calc_bin(dih[index][j], n, core_frac);
532 accum = n * accum + b - 1;
547 chi_prtrj[j] = accum;
559 /* print cuml rotamer vs time */
560 print_one(oenv, "chiproduct", dlist[i].name, "chi product for",
561 "cumulative rotamer", nframes, time, chi_prtrj);
564 /* make a histogram pf culm. rotamer occupancy too */
565 snew(chi_prhist, nbin);
566 make_histo(NULL, nframes, chi_prtrj, nbin, chi_prhist, 0, nbin);
569 sprintf(hisfile, "histo-chiprod%s.xvg", dlist[i].name);
570 sprintf(histitle, "cumulative rotamer distribution for %s", dlist[i].name);
571 fprintf(stderr, " and %s ", hisfile);
572 fp = xvgropen(hisfile, histitle, "number", "", oenv);
573 if (output_env_get_print_xvgr_codes(oenv))
575 fprintf(fp, "@ xaxis tick on\n");
576 fprintf(fp, "@ xaxis tick major 1\n");
577 fprintf(fp, "@ type xy\n");
579 for (k = 0; (k < nbin); k++)
583 fprintf(fp, "%5d %10g\n", k, (1.0*chi_prhist[k])/nframes);
587 fprintf(fp, "%5d %10d\n", k, chi_prhist[k]);
590 fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
594 /* and finally print out occupancies to a single file */
595 /* get the gmx from-1 res nr by setting a ptr to the number part
596 * of dlist[i].name - potential bug for 4-letter res names... */
597 namept = dlist[i].name + 3;
598 fprintf(fpall, "%5s ", namept);
599 for (k = 0; (k < nbin); k++)
603 fprintf(fpall, " %10g", (1.0*chi_prhist[k])/nframes);
607 fprintf(fpall, " %10d", chi_prhist[k]);
610 fprintf(fpall, "\n");
619 fprintf(stderr, "\n");
623 void calc_distribution_props(int nh, int histo[], real start,
624 int nkkk, t_karplus kkk[],
627 real d, dc, ds, c1, c2, tdc, tds;
628 real fac, ang, invth, Jc;
633 gmx_fatal(FARGS, "No points in histogram (%s, %d)", __FILE__, __LINE__);
637 /* Compute normalisation factor */
639 for (j = 0; (j < nh); j++)
645 for (i = 0; (i < nkkk); i++)
651 for (j = 0; (j < nh); j++)
661 for (i = 0; (i < nkkk); i++)
663 c1 = cos(ang+kkk[i].offset);
665 Jc = (kkk[i].A*c2 + kkk[i].B*c1 + kkk[i].C);
666 kkk[i].Jc += histo[j]*Jc;
667 kkk[i].Jcsig += histo[j]*sqr(Jc);
670 for (i = 0; (i < nkkk); i++)
673 kkk[i].Jcsig = sqrt(kkk[i].Jcsig/th-sqr(kkk[i].Jc));
675 *S2 = tdc*tdc+tds*tds;
678 static void calc_angles(struct t_pbc *pbc,
679 int n3, atom_id index[], real ang[], rvec x_s[])
685 for (i = ix = 0; (ix < n3); i++, ix += 3)
687 ang[i] = bond_angle(x_s[index[ix]], x_s[index[ix+1]], x_s[index[ix+2]],
688 pbc, r_ij, r_kj, &costh, &t1, &t2);
692 fprintf(debug, "Angle[0]=%g, costh=%g, index0 = %d, %d, %d\n",
693 ang[0], costh, index[0], index[1], index[2]);
694 pr_rvec(debug, 0, "rij", r_ij, DIM, TRUE);
695 pr_rvec(debug, 0, "rkj", r_kj, DIM, TRUE);
699 static real calc_fraction(real angles[], int nangles)
702 real trans = 0, gauche = 0;
705 for (i = 0; i < nangles; i++)
707 angle = angles[i] * RAD2DEG;
709 if (angle > 135 && angle < 225)
713 else if (angle > 270 && angle < 330)
717 else if (angle < 90 && angle > 30)
722 if (trans+gauche > 0)
724 return trans/(trans+gauche);
732 static void calc_dihs(struct t_pbc *pbc,
733 int n4, atom_id index[], real ang[], rvec x_s[])
735 int i, ix, t1, t2, t3;
736 rvec r_ij, r_kj, r_kl, m, n;
739 for (i = ix = 0; (ix < n4); i++, ix += 4)
741 aaa = dih_angle(x_s[index[ix]], x_s[index[ix+1]], x_s[index[ix+2]],
742 x_s[index[ix+3]], pbc,
743 r_ij, r_kj, r_kl, m, n,
744 &sign, &t1, &t2, &t3);
746 ang[i] = aaa; /* not taking into account ryckaert bellemans yet */
750 void make_histo(FILE *log,
751 int ndata, real data[], int npoints, int histo[],
752 real minx, real maxx)
759 minx = maxx = data[0];
760 for (i = 1; (i < ndata); i++)
762 minx = min(minx, data[i]);
763 maxx = max(maxx, data[i]);
765 fprintf(log, "Min data: %10g Max data: %10g\n", minx, maxx);
767 dx = (double)npoints/(maxx-minx);
771 "Histogramming: ndata=%d, nhisto=%d, minx=%g,maxx=%g,dx=%g\n",
772 ndata, npoints, minx, maxx, dx);
774 for (i = 0; (i < ndata); i++)
776 ind = (data[i]-minx)*dx;
777 if ((ind >= 0) && (ind < npoints))
783 fprintf(log, "index = %d, data[%d] = %g\n", ind, i, data[i]);
788 void normalize_histo(int npoints, int histo[], real dx, real normhisto[])
794 for (i = 0; (i < npoints); i++)
800 fprintf(stderr, "Empty histogram!\n");
804 for (i = 0; (i < npoints); i++)
806 normhisto[i] = fac*histo[i];
810 void read_ang_dih(const char *trj_fn,
811 gmx_bool bAngles, gmx_bool bSaveAll, gmx_bool bRb, gmx_bool bPBC,
812 int maxangstat, int angstat[],
813 int *nframes, real **time,
814 int isize, atom_id index[],
818 const output_env_t oenv)
822 int i, angind, natoms, total, teller;
823 int nangles, n_alloc;
824 real t, fraction, pifac, aa, angle;
832 natoms = read_first_x(oenv, &status, trj_fn, &t, &x, box);
844 snew(angles[cur], nangles);
845 snew(angles[prev], nangles);
847 /* Start the loop over frames */
857 if (teller >= n_alloc)
862 for (i = 0; (i < nangles); i++)
864 srenew(dih[i], n_alloc);
867 srenew(*time, n_alloc);
868 srenew(*trans_frac, n_alloc);
869 srenew(*aver_angle, n_alloc);
876 set_pbc(pbc, -1, box);
881 calc_angles(pbc, isize, index, angles[cur], x);
885 calc_dihs(pbc, isize, index, angles[cur], x);
888 fraction = calc_fraction(angles[cur], nangles);
889 (*trans_frac)[teller] = fraction;
891 /* Change Ryckaert-Bellemans dihedrals to polymer convention
892 * Modified 990913 by Erik:
893 * We actually shouldn't change the convention, since it's
894 * calculated from polymer above, but we change the intervall
895 * from [-180,180] to [0,360].
899 for (i = 0; (i < nangles); i++)
901 if (angles[cur][i] <= 0.0)
903 angles[cur][i] += 2*M_PI;
908 /* Periodicity in dihedral space... */
911 for (i = 0; (i < nangles); i++)
913 real dd = angles[cur][i];
914 angles[cur][i] = atan2(sin(dd), cos(dd));
921 for (i = 0; (i < nangles); i++)
923 while (angles[cur][i] <= angles[prev][i] - M_PI)
925 angles[cur][i] += 2*M_PI;
927 while (angles[cur][i] > angles[prev][i] + M_PI)
929 angles[cur][i] -= 2*M_PI;
938 for (i = 0; (i < nangles); i++)
940 aa = aa+angles[cur][i];
942 /* angle in rad / 2Pi * max determines bin. bins go from 0 to maxangstat,
943 even though scale goes from -pi to pi (dihedral) or -pi/2 to pi/2
944 (angle) Basically: translate the x-axis by Pi. Translate it back by
948 angle = angles[cur][i];
951 while (angle < -M_PI)
955 while (angle >= M_PI)
963 /* Update the distribution histogram */
964 angind = (int) ((angle*maxangstat)/pifac + 0.5);
965 if (angind == maxangstat)
969 if ( (angind < 0) || (angind >= maxangstat) )
971 /* this will never happen */
972 gmx_fatal(FARGS, "angle (%f) index out of range (0..%d) : %d\n",
973 angle, maxangstat, angind);
977 if (angind == maxangstat)
979 fprintf(stderr, "angle %d fr %d = %g\n", i, cur, angle);
985 /* average over all angles */
986 (*aver_angle)[teller] = (aa/nangles);
988 /* this copies all current dih. angles to dih[i], teller is frame */
991 for (i = 0; i < nangles; i++)
993 dih[i][teller] = angles[cur][i];
1000 /* Increment loop counter */
1003 while (read_next_x(oenv, status, &t, x, box));
1008 sfree(angles[prev]);