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.
45 #include "gromacs/utility/smalloc.h"
53 #include "gromacs/fileio/confio.h"
54 #include "gromacs/fileio/trxio.h"
56 #include "gromacs/legacyheaders/gmx_fatal.h"
58 void print_one(const output_env_t oenv, const char *base, const char *name,
59 const char *title, const char *ylabel, int nf, real time[],
63 char buf[256], t2[256];
66 sprintf(buf, "%s%s.xvg", base, name);
67 fprintf(stderr, "\rPrinting %s ", buf);
68 sprintf(t2, "%s %s", title, name);
69 fp = xvgropen(buf, t2, "Time (ps)", ylabel, oenv);
70 for (k = 0; (k < nf); k++)
72 fprintf(fp, "%10g %10g\n", time[k], data[k]);
77 static int calc_RBbin(real phi, int gmx_unused multiplicity, real gmx_unused core_frac)
79 /* multiplicity and core_frac NOT used,
80 * just given to enable use of pt-to-fn in caller low_ana_dih_trans*/
81 static const real r30 = M_PI/6.0;
82 static const real r90 = M_PI/2.0;
83 static const real r150 = M_PI*5.0/6.0;
85 if ((phi < r30) && (phi > -r30))
89 else if ((phi > -r150) && (phi < -r90))
93 else if ((phi < r150) && (phi > r90))
100 static int calc_Nbin(real phi, int multiplicity, real core_frac)
102 static const real r360 = 360*DEG2RAD;
103 real rot_width, core_width, core_offset, low, hi;
105 /* with multiplicity 3 and core_frac 0.5
106 * 0<g(-)<120, 120<t<240, 240<g(+)<360
107 * 0< bin0 < 30, 30<bin1<90, 90<bin0<150, 150<bin2<210, 210<bin0<270, 270<bin3<330, 330<bin0<360
108 * so with multiplicity 3, bin1 is core g(-), bin2 is core t, bin3 is
109 core g(+), bin0 is between rotamers */
115 rot_width = 360/multiplicity;
116 core_width = core_frac * rot_width;
117 core_offset = (rot_width - core_width)/2.0;
118 for (bin = 1; bin <= multiplicity; bin++)
120 low = ((bin - 1) * rot_width ) + core_offset;
121 hi = ((bin - 1) * rot_width ) + core_offset + core_width;
124 if ((phi > low) && (phi < hi))
132 void ana_dih_trans(const char *fn_trans, const char *fn_histo,
133 real **dih, int nframes, int nangles,
134 const char *grpname, real *time, gmx_bool bRb,
135 const output_env_t oenv)
137 /* just a wrapper; declare extra args, then chuck away at end. */
145 snew(multiplicity, nangles);
146 for (k = 0; (k < nangles); k++)
151 low_ana_dih_trans(TRUE, fn_trans, TRUE, fn_histo, maxchi,
152 dih, nlist, dlist, nframes,
153 nangles, grpname, multiplicity, time, bRb, 0.5, oenv);
159 void low_ana_dih_trans(gmx_bool bTrans, const char *fn_trans,
160 gmx_bool bHisto, const char *fn_histo, int maxchi,
161 real **dih, int nlist, t_dlist dlist[], int nframes,
162 int nangles, const char *grpname, int multiplicity[],
163 real *time, gmx_bool bRb, real core_frac,
164 const output_env_t oenv)
169 int i, j, k, Dih, ntrans;
170 int cur_bin, new_bin;
173 int (*calc_bin)(real, int, real);
180 /* Assumes the frames are equally spaced in time */
181 dt = (time[nframes-1]-time[0])/(nframes-1);
183 /* Analysis of dihedral transitions */
184 fprintf(stderr, "Now calculating transitions...\n");
188 calc_bin = calc_RBbin;
192 calc_bin = calc_Nbin;
195 for (k = 0; k < NROT; k++)
197 snew(rot_occ[k], nangles);
198 for (i = 0; (i < nangles); i++)
206 /* dih[i][j] is the dihedral angle i in frame j */
208 for (i = 0; (i < nangles); i++)
213 mind = maxd = prev = dih[i][0];
215 cur_bin = calc_bin(dih[i][0], multiplicity[i], core_frac);
216 rot_occ[cur_bin][i]++;
218 for (j = 1; (j < nframes); j++)
220 new_bin = calc_bin(dih[i][j], multiplicity[i], core_frac);
221 rot_occ[new_bin][i]++;
227 else if ((new_bin != 0) && (cur_bin != new_bin))
235 /* why is all this md rubbish periodic? Remove 360 degree periodicity */
236 if ( (dih[i][j] - prev) > M_PI)
240 else if ( (dih[i][j] - prev) < -M_PI)
247 mind = min(mind, dih[i][j]);
248 maxd = max(maxd, dih[i][j]);
249 if ( (maxd - mind) > 2*M_PI/3) /* or 120 degrees, assuming */
250 { /* multiplicity 3. Not so general.*/
253 maxd = mind = dih[i][j]; /* get ready for next transition */
258 for (k = 0; k < NROT; k++)
260 rot_occ[k][i] /= nframes;
263 fprintf(stderr, "Total number of transitions: %10d\n", ntrans);
266 ttime = (dt*nframes*nangles)/ntrans;
267 fprintf(stderr, "Time between transitions: %10.3f ps\n", ttime);
270 /* new by grs - copy transitions from tr_h[] to dlist->ntr[]
271 * and rotamer populations from rot_occ to dlist->rot_occ[]
272 * based on fn histogramming in g_chi. diff roles for i and j here */
275 for (Dih = 0; (Dih < NONCHI+maxchi); Dih++)
277 for (i = 0; (i < nlist); i++)
279 if (((Dih < edOmega) ) ||
280 ((Dih == edOmega) && (has_dihedral(edOmega, &(dlist[i])))) ||
281 ((Dih > edOmega) && (dlist[i].atm.Cn[Dih-NONCHI+3] != -1)))
283 /* grs debug printf("Not OK? i %d j %d Dih %d \n", i, j, Dih) ; */
284 dlist[i].ntr[Dih] = tr_h[j];
285 for (k = 0; k < NROT; k++)
287 dlist[i].rot_occ[Dih][k] = rot_occ[k][j];
294 /* end addition by grs */
298 sprintf(title, "Number of transitions: %s", grpname);
299 fp = xvgropen(fn_trans, title, "Time (ps)", "# transitions/timeframe", oenv);
300 for (j = 0; (j < nframes); j++)
302 fprintf(fp, "%10.3f %10d\n", time[j], tr_f[j]);
307 /* Compute histogram from # transitions per dihedral */
309 for (j = 0; (j < nframes); j++)
313 for (i = 0; (i < nangles); i++)
317 for (j = nframes; ((tr_f[j-1] == 0) && (j > 0)); j--)
325 sprintf(title, "Transition time: %s", grpname);
326 fp = xvgropen(fn_histo, title, "Time (ps)", "#", oenv);
327 for (i = j-1; (i > 0); i--)
331 fprintf(fp, "%10.3f %10d\n", ttime/i, tr_f[i]);
339 for (k = 0; k < NROT; k++)
346 void mk_multiplicity_lookup (int *multiplicity, int maxchi,
347 int nlist, t_dlist dlist[], int nangles)
349 /* new by grs - for dihedral j (as in dih[j]) get multiplicity from dlist
350 * and store in multiplicity[j]
357 for (Dih = 0; (Dih < NONCHI+maxchi); Dih++)
359 for (i = 0; (i < nlist); i++)
361 strncpy(name, dlist[i].name, 3);
363 if (((Dih < edOmega) ) ||
364 ((Dih == edOmega) && (has_dihedral(edOmega, &(dlist[i])))) ||
365 ((Dih > edOmega) && (dlist[i].atm.Cn[Dih-NONCHI+3] != -1)))
367 /* default - we will correct the rest below */
370 /* make omegas 2fold, though doesn't make much more sense than 3 */
371 if (Dih == edOmega && (has_dihedral(edOmega, &(dlist[i]))))
376 /* dihedrals to aromatic rings, COO, CONH2 or guanidinium are 2fold*/
377 if (Dih > edOmega && (dlist[i].atm.Cn[Dih-NONCHI+3] != -1))
379 if ( ((strstr(name, "PHE") != NULL) && (Dih == edChi2)) ||
380 ((strstr(name, "TYR") != NULL) && (Dih == edChi2)) ||
381 ((strstr(name, "PTR") != NULL) && (Dih == edChi2)) ||
382 ((strstr(name, "TRP") != NULL) && (Dih == edChi2)) ||
383 ((strstr(name, "HIS") != NULL) && (Dih == edChi2)) ||
384 ((strstr(name, "GLU") != NULL) && (Dih == edChi3)) ||
385 ((strstr(name, "ASP") != NULL) && (Dih == edChi2)) ||
386 ((strstr(name, "GLN") != NULL) && (Dih == edChi3)) ||
387 ((strstr(name, "ASN") != NULL) && (Dih == edChi2)) ||
388 ((strstr(name, "ARG") != NULL) && (Dih == edChi4)) )
399 fprintf(stderr, "WARNING: not all dihedrals found in topology (only %d out of %d)!\n",
402 /* Check for remaining dihedrals */
403 for (; (j < nangles); j++)
410 void mk_chi_lookup (int **lookup, int maxchi,
411 int nlist, t_dlist dlist[])
414 /* by grs. should rewrite everything to use this. (but haven't,
415 * and at mmt only used in get_chi_product_traj
416 * returns the dihed number given the residue number (from-0)
417 * and chi (from-0) nr. -1 for chi undefined for that res (eg gly, ala..)*/
422 /* NONCHI points to chi1, therefore we have to start counting there. */
423 for (Dih = NONCHI; (Dih < NONCHI+maxchi); Dih++)
425 for (i = 0; (i < nlist); i++)
428 if (((Dih < edOmega) ) ||
429 ((Dih == edOmega) && (has_dihedral(edOmega, &(dlist[i])))) ||
430 ((Dih > edOmega) && (dlist[i].atm.Cn[Dih-NONCHI+3] != -1)))
432 /* grs debug printf("Not OK? i %d j %d Dih %d \n", i, j, Dih) ; */
449 void get_chi_product_traj (real **dih, int nframes, int nlist,
450 int maxchi, t_dlist dlist[], real time[],
451 int **lookup, int *multiplicity, gmx_bool bRb, gmx_bool bNormalize,
452 real core_frac, gmx_bool bAll, const char *fnall,
453 const output_env_t oenv)
456 gmx_bool bRotZero, bHaveChi = FALSE;
457 int accum = 0, index, i, j, k, Xi, n, b;
462 char hisfile[256], histitle[256], *namept;
464 int (*calc_bin)(real, int, real);
466 /* Analysis of dihedral transitions */
467 fprintf(stderr, "Now calculating Chi product trajectories...\n");
471 calc_bin = calc_RBbin;
475 calc_bin = calc_Nbin;
478 snew(chi_prtrj, nframes);
480 /* file for info on all residues */
483 fpall = xvgropen(fnall, "Cumulative Rotamers", "Residue", "Probability", oenv);
487 fpall = xvgropen(fnall, "Cumulative Rotamers", "Residue", "# Counts", oenv);
490 for (i = 0; (i < nlist); i++)
493 /* get nbin, the nr. of cumulative rotamers that need to be considered */
495 for (Xi = 0; Xi < maxchi; Xi++)
497 index = lookup[i][Xi]; /* chi_(Xi+1) of res i (-1 if off end) */
500 n = multiplicity[index];
504 nbin += 1; /* for the "zero rotamer", outside the core region */
506 for (j = 0; (j < nframes); j++)
511 index = lookup[i][0]; /* index into dih of chi1 of res i */
520 b = calc_bin(dih[index][j], multiplicity[index], core_frac);
526 for (Xi = 1; Xi < maxchi; Xi++)
528 index = lookup[i][Xi]; /* chi_(Xi+1) of res i (-1 if off end) */
531 n = multiplicity[index];
532 b = calc_bin(dih[index][j], n, core_frac);
533 accum = n * accum + b - 1;
548 chi_prtrj[j] = accum;
560 /* print cuml rotamer vs time */
561 print_one(oenv, "chiproduct", dlist[i].name, "chi product for",
562 "cumulative rotamer", nframes, time, chi_prtrj);
565 /* make a histogram pf culm. rotamer occupancy too */
566 snew(chi_prhist, nbin);
567 make_histo(NULL, nframes, chi_prtrj, nbin, chi_prhist, 0, nbin);
570 sprintf(hisfile, "histo-chiprod%s.xvg", dlist[i].name);
571 sprintf(histitle, "cumulative rotamer distribution for %s", dlist[i].name);
572 fprintf(stderr, " and %s ", hisfile);
573 fp = xvgropen(hisfile, histitle, "number", "", oenv);
574 if (output_env_get_print_xvgr_codes(oenv))
576 fprintf(fp, "@ xaxis tick on\n");
577 fprintf(fp, "@ xaxis tick major 1\n");
578 fprintf(fp, "@ type xy\n");
580 for (k = 0; (k < nbin); k++)
584 fprintf(fp, "%5d %10g\n", k, (1.0*chi_prhist[k])/nframes);
588 fprintf(fp, "%5d %10d\n", k, chi_prhist[k]);
591 fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
595 /* and finally print out occupancies to a single file */
596 /* get the gmx from-1 res nr by setting a ptr to the number part
597 * of dlist[i].name - potential bug for 4-letter res names... */
598 namept = dlist[i].name + 3;
599 fprintf(fpall, "%5s ", namept);
600 for (k = 0; (k < nbin); k++)
604 fprintf(fpall, " %10g", (1.0*chi_prhist[k])/nframes);
608 fprintf(fpall, " %10d", chi_prhist[k]);
611 fprintf(fpall, "\n");
620 fprintf(stderr, "\n");
624 void calc_distribution_props(int nh, int histo[], real start,
625 int nkkk, t_karplus kkk[],
628 real d, dc, ds, c1, c2, tdc, tds;
629 real fac, ang, invth, Jc;
634 gmx_fatal(FARGS, "No points in histogram (%s, %d)", __FILE__, __LINE__);
638 /* Compute normalisation factor */
640 for (j = 0; (j < nh); j++)
646 for (i = 0; (i < nkkk); i++)
652 for (j = 0; (j < nh); j++)
662 for (i = 0; (i < nkkk); i++)
664 c1 = cos(ang+kkk[i].offset);
666 Jc = (kkk[i].A*c2 + kkk[i].B*c1 + kkk[i].C);
667 kkk[i].Jc += histo[j]*Jc;
668 kkk[i].Jcsig += histo[j]*sqr(Jc);
671 for (i = 0; (i < nkkk); i++)
674 kkk[i].Jcsig = sqrt(kkk[i].Jcsig/th-sqr(kkk[i].Jc));
676 *S2 = tdc*tdc+tds*tds;
679 static void calc_angles(t_pbc *pbc,
680 int n3, atom_id index[], real ang[], rvec x_s[])
686 for (i = ix = 0; (ix < n3); i++, ix += 3)
688 ang[i] = bond_angle(x_s[index[ix]], x_s[index[ix+1]], x_s[index[ix+2]],
689 pbc, r_ij, r_kj, &costh, &t1, &t2);
693 fprintf(debug, "Angle[0]=%g, costh=%g, index0 = %d, %d, %d\n",
694 ang[0], costh, index[0], index[1], index[2]);
695 pr_rvec(debug, 0, "rij", r_ij, DIM, TRUE);
696 pr_rvec(debug, 0, "rkj", r_kj, DIM, TRUE);
700 static real calc_fraction(real angles[], int nangles)
703 real trans = 0, gauche = 0;
706 for (i = 0; i < nangles; i++)
708 angle = angles[i] * RAD2DEG;
710 if (angle > 135 && angle < 225)
714 else if (angle > 270 && angle < 330)
718 else if (angle < 90 && angle > 30)
723 if (trans+gauche > 0)
725 return trans/(trans+gauche);
733 static void calc_dihs(t_pbc *pbc,
734 int n4, atom_id index[], real ang[], rvec x_s[])
736 int i, ix, t1, t2, t3;
737 rvec r_ij, r_kj, r_kl, m, n;
740 for (i = ix = 0; (ix < n4); i++, ix += 4)
742 aaa = dih_angle(x_s[index[ix]], x_s[index[ix+1]], x_s[index[ix+2]],
743 x_s[index[ix+3]], pbc,
744 r_ij, r_kj, r_kl, m, n,
745 &sign, &t1, &t2, &t3);
747 ang[i] = aaa; /* not taking into account ryckaert bellemans yet */
751 void make_histo(FILE *log,
752 int ndata, real data[], int npoints, int histo[],
753 real minx, real maxx)
760 minx = maxx = data[0];
761 for (i = 1; (i < ndata); i++)
763 minx = min(minx, data[i]);
764 maxx = max(maxx, data[i]);
766 fprintf(log, "Min data: %10g Max data: %10g\n", minx, maxx);
768 dx = (double)npoints/(maxx-minx);
772 "Histogramming: ndata=%d, nhisto=%d, minx=%g,maxx=%g,dx=%g\n",
773 ndata, npoints, minx, maxx, dx);
775 for (i = 0; (i < ndata); i++)
777 ind = (data[i]-minx)*dx;
778 if ((ind >= 0) && (ind < npoints))
784 fprintf(log, "index = %d, data[%d] = %g\n", ind, i, data[i]);
789 void normalize_histo(int npoints, int histo[], real dx, real normhisto[])
795 for (i = 0; (i < npoints); i++)
801 fprintf(stderr, "Empty histogram!\n");
805 for (i = 0; (i < npoints); i++)
807 normhisto[i] = fac*histo[i];
811 void read_ang_dih(const char *trj_fn,
812 gmx_bool bAngles, gmx_bool bSaveAll, gmx_bool bRb, gmx_bool bPBC,
813 int maxangstat, int angstat[],
814 int *nframes, real **time,
815 int isize, atom_id index[],
819 const output_env_t oenv)
823 int i, angind, natoms, total, teller;
824 int nangles, n_alloc;
825 real t, fraction, pifac, aa, angle;
833 natoms = read_first_x(oenv, &status, trj_fn, &t, &x, box);
845 snew(angles[cur], nangles);
846 snew(angles[prev], nangles);
848 /* Start the loop over frames */
858 if (teller >= n_alloc)
863 for (i = 0; (i < nangles); i++)
865 srenew(dih[i], n_alloc);
868 srenew(*time, n_alloc);
869 srenew(*trans_frac, n_alloc);
870 srenew(*aver_angle, n_alloc);
877 set_pbc(pbc, -1, box);
882 calc_angles(pbc, isize, index, angles[cur], x);
886 calc_dihs(pbc, isize, index, angles[cur], x);
889 fraction = calc_fraction(angles[cur], nangles);
890 (*trans_frac)[teller] = fraction;
892 /* Change Ryckaert-Bellemans dihedrals to polymer convention
893 * Modified 990913 by Erik:
894 * We actually shouldn't change the convention, since it's
895 * calculated from polymer above, but we change the intervall
896 * from [-180,180] to [0,360].
900 for (i = 0; (i < nangles); i++)
902 if (angles[cur][i] <= 0.0)
904 angles[cur][i] += 2*M_PI;
909 /* Periodicity in dihedral space... */
912 for (i = 0; (i < nangles); i++)
914 real dd = angles[cur][i];
915 angles[cur][i] = atan2(sin(dd), cos(dd));
922 for (i = 0; (i < nangles); i++)
924 while (angles[cur][i] <= angles[prev][i] - M_PI)
926 angles[cur][i] += 2*M_PI;
928 while (angles[cur][i] > angles[prev][i] + M_PI)
930 angles[cur][i] -= 2*M_PI;
939 for (i = 0; (i < nangles); i++)
941 aa = aa+angles[cur][i];
943 /* angle in rad / 2Pi * max determines bin. bins go from 0 to maxangstat,
944 even though scale goes from -pi to pi (dihedral) or -pi/2 to pi/2
945 (angle) Basically: translate the x-axis by Pi. Translate it back by
949 angle = angles[cur][i];
952 while (angle < -M_PI)
956 while (angle >= M_PI)
964 /* Update the distribution histogram */
965 angind = (int) ((angle*maxangstat)/pifac + 0.5);
966 if (angind == maxangstat)
970 if ( (angind < 0) || (angind >= maxangstat) )
972 /* this will never happen */
973 gmx_fatal(FARGS, "angle (%f) index out of range (0..%d) : %d\n",
974 angle, maxangstat, angind);
978 if (angind == maxangstat)
980 fprintf(stderr, "angle %d fr %d = %g\n", i, cur, angle);
986 /* average over all angles */
987 (*aver_angle)[teller] = (aa/nangles);
989 /* this copies all current dih. angles to dih[i], teller is frame */
992 for (i = 0; i < nangles; i++)
994 dih[i][teller] = angles[cur][i];
1001 /* Increment loop counter */
1004 while (read_next_x(oenv, status, &t, x, box));
1009 sfree(angles[prev]);