Cleaned up memory usage in gmx traj and trjconv
[gromacs.git] / src / gromacs / gmxana / anadih.cpp
blob2393f556f012acf64784ff3e536f7cbd27c3ee85
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, 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.
37 #include "gmxpre.h"
39 #include <cmath>
40 #include <cstdio>
41 #include <cstring>
43 #include <algorithm>
45 #include "gromacs/fileio/confio.h"
46 #include "gromacs/fileio/trxio.h"
47 #include "gromacs/fileio/xvgr.h"
48 #include "gromacs/gmxana/gstat.h"
49 #include "gromacs/listed-forces/bonded.h"
50 #include "gromacs/math/functions.h"
51 #include "gromacs/math/units.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/math/vecdump.h"
54 #include "gromacs/pbcutil/pbc.h"
55 #include "gromacs/utility/fatalerror.h"
56 #include "gromacs/utility/smalloc.h"
58 void print_one(const gmx_output_env_t *oenv, const char *base, const char *name,
59 const char *title, const char *ylabel, int nf, real time[],
60 real data[])
62 FILE *fp;
63 char buf[256], t2[256];
64 int k;
66 sprintf(buf, "%s%s.xvg", base, name);
67 fprintf(stderr, "\rPrinting %s ", buf);
68 fflush(stderr);
69 sprintf(t2, "%s %s", title, name);
70 fp = xvgropen(buf, t2, "Time (ps)", ylabel, oenv);
71 for (k = 0; (k < nf); k++)
73 fprintf(fp, "%10g %10g\n", time[k], data[k]);
75 xvgrclose(fp);
78 static int calc_RBbin(real phi, int gmx_unused multiplicity, real gmx_unused core_frac)
80 /* multiplicity and core_frac NOT used,
81 * just given to enable use of pt-to-fn in caller low_ana_dih_trans*/
82 static const real r30 = M_PI/6.0;
83 static const real r90 = M_PI/2.0;
84 static const real r150 = M_PI*5.0/6.0;
86 if ((phi < r30) && (phi > -r30))
88 return 1;
90 else if ((phi > -r150) && (phi < -r90))
92 return 2;
94 else if ((phi < r150) && (phi > r90))
96 return 3;
98 return 0;
101 static int calc_Nbin(real phi, int multiplicity, real core_frac)
103 static const real r360 = 360*DEG2RAD;
104 real rot_width, core_width, core_offset, low, hi;
105 int bin;
106 /* with multiplicity 3 and core_frac 0.5
107 * 0<g(-)<120, 120<t<240, 240<g(+)<360
108 * 0< bin0 < 30, 30<bin1<90, 90<bin0<150, 150<bin2<210, 210<bin0<270, 270<bin3<330, 330<bin0<360
109 * so with multiplicity 3, bin1 is core g(-), bin2 is core t, bin3 is
110 core g(+), bin0 is between rotamers */
111 if (phi < 0)
113 phi += r360;
116 rot_width = 360/multiplicity;
117 core_width = core_frac * rot_width;
118 core_offset = (rot_width - core_width)/2.0;
119 for (bin = 1; bin <= multiplicity; bin++)
121 low = ((bin - 1) * rot_width ) + core_offset;
122 hi = ((bin - 1) * rot_width ) + core_offset + core_width;
123 low *= DEG2RAD;
124 hi *= DEG2RAD;
125 if ((phi > low) && (phi < hi))
127 return bin;
130 return 0;
133 void ana_dih_trans(const char *fn_trans, const char *fn_histo,
134 real **dih, int nframes, int nangles,
135 const char *grpname, real *time, gmx_bool bRb,
136 const gmx_output_env_t *oenv)
138 /* just a wrapper; declare extra args, then chuck away at end. */
139 int maxchi = 0;
140 t_dlist *dlist;
141 int *multiplicity;
142 int nlist = nangles;
143 int k;
145 snew(dlist, nlist);
146 snew(multiplicity, nangles);
147 for (k = 0; (k < nangles); k++)
149 multiplicity[k] = 3;
152 low_ana_dih_trans(TRUE, fn_trans, TRUE, fn_histo, maxchi,
153 dih, nlist, dlist, nframes,
154 nangles, grpname, multiplicity, time, bRb, 0.5, oenv);
155 sfree(dlist);
156 sfree(multiplicity);
160 void low_ana_dih_trans(gmx_bool bTrans, const char *fn_trans,
161 gmx_bool bHisto, const char *fn_histo, int maxchi,
162 real **dih, int nlist, t_dlist dlist[], int nframes,
163 int nangles, const char *grpname, int multiplicity[],
164 real *time, gmx_bool bRb, real core_frac,
165 const gmx_output_env_t *oenv)
167 FILE *fp;
168 int *tr_f, *tr_h;
169 char title[256];
170 int i, j, k, Dih, ntrans;
171 int cur_bin, new_bin;
172 real ttime;
173 real *rot_occ[NROT];
174 int (*calc_bin)(real, int, real);
175 real dt;
177 if (1 <= nframes)
179 return;
181 /* Assumes the frames are equally spaced in time */
182 dt = (time[nframes-1]-time[0])/(nframes-1);
184 /* Analysis of dihedral transitions */
185 fprintf(stderr, "Now calculating transitions...\n");
187 if (bRb)
189 calc_bin = calc_RBbin;
191 else
193 calc_bin = calc_Nbin;
196 for (k = 0; k < NROT; k++)
198 snew(rot_occ[k], nangles);
199 for (i = 0; (i < nangles); i++)
201 rot_occ[k][i] = 0;
204 snew(tr_h, nangles);
205 snew(tr_f, nframes);
207 /* dih[i][j] is the dihedral angle i in frame j */
208 ntrans = 0;
209 for (i = 0; (i < nangles); i++)
212 /*#define OLDIE*/
213 #ifdef OLDIE
214 mind = maxd = prev = dih[i][0];
215 #else
216 cur_bin = calc_bin(dih[i][0], multiplicity[i], core_frac);
217 rot_occ[cur_bin][i]++;
218 #endif
219 for (j = 1; (j < nframes); j++)
221 new_bin = calc_bin(dih[i][j], multiplicity[i], core_frac);
222 rot_occ[new_bin][i]++;
223 #ifndef OLDIE
224 if (cur_bin == 0)
226 cur_bin = new_bin;
228 else if ((new_bin != 0) && (cur_bin != new_bin))
230 cur_bin = new_bin;
231 tr_f[j]++;
232 tr_h[i]++;
233 ntrans++;
235 #else
236 /* why is all this md rubbish periodic? Remove 360 degree periodicity */
237 if ( (dih[i][j] - prev) > M_PI)
239 dih[i][j] -= 2*M_PI;
241 else if ( (dih[i][j] - prev) < -M_PI)
243 dih[i][j] += 2*M_PI;
246 prev = dih[i][j];
248 mind = std::min(mind, dih[i][j]);
249 maxd = std::max(maxd, dih[i][j]);
250 if ( (maxd - mind) > 2*M_PI/3) /* or 120 degrees, assuming */
251 { /* multiplicity 3. Not so general.*/
252 tr_f[j]++;
253 tr_h[i]++;
254 maxd = mind = dih[i][j]; /* get ready for next transition */
255 ntrans++;
257 #endif
258 } /* end j */
259 for (k = 0; k < NROT; k++)
261 rot_occ[k][i] /= nframes;
263 } /* end i */
264 fprintf(stderr, "Total number of transitions: %10d\n", ntrans);
265 if (ntrans > 0)
267 ttime = (dt*nframes*nangles)/ntrans;
268 fprintf(stderr, "Time between transitions: %10.3f ps\n", ttime);
271 /* new by grs - copy transitions from tr_h[] to dlist->ntr[]
272 * and rotamer populations from rot_occ to dlist->rot_occ[]
273 * based on fn histogramming in g_chi. diff roles for i and j here */
275 j = 0;
276 for (Dih = 0; (Dih < NONCHI+maxchi); Dih++)
278 for (i = 0; (i < nlist); i++)
280 if (((Dih < edOmega) ) ||
281 ((Dih == edOmega) && (has_dihedral(edOmega, &(dlist[i])))) ||
282 ((Dih > edOmega) && (dlist[i].atm.Cn[Dih-NONCHI+3] != -1)))
284 /* grs debug printf("Not OK? i %d j %d Dih %d \n", i, j, Dih) ; */
285 dlist[i].ntr[Dih] = tr_h[j];
286 for (k = 0; k < NROT; k++)
288 dlist[i].rot_occ[Dih][k] = rot_occ[k][j];
290 j++;
295 /* end addition by grs */
297 if (bTrans)
299 sprintf(title, "Number of transitions: %s", grpname);
300 fp = xvgropen(fn_trans, title, "Time (ps)", "# transitions/timeframe", oenv);
301 for (j = 0; (j < nframes); j++)
303 fprintf(fp, "%10.3f %10d\n", time[j], tr_f[j]);
305 xvgrclose(fp);
308 /* Compute histogram from # transitions per dihedral */
309 /* Use old array */
310 for (j = 0; (j < nframes); j++)
312 tr_f[j] = 0;
314 for (i = 0; (i < nangles); i++)
316 tr_f[tr_h[i]]++;
318 for (j = nframes; ((tr_f[j-1] == 0) && (j > 0)); j--)
323 ttime = dt*nframes;
324 if (bHisto)
326 sprintf(title, "Transition time: %s", grpname);
327 fp = xvgropen(fn_histo, title, "Time (ps)", "#", oenv);
328 for (i = j-1; (i > 0); i--)
330 if (tr_f[i] != 0)
332 fprintf(fp, "%10.3f %10d\n", ttime/i, tr_f[i]);
335 xvgrclose(fp);
338 sfree(tr_f);
339 sfree(tr_h);
340 for (k = 0; k < NROT; k++)
342 sfree(rot_occ[k]);
347 void mk_multiplicity_lookup (int *multiplicity, int maxchi,
348 int nlist, t_dlist dlist[], int nangles)
350 /* new by grs - for dihedral j (as in dih[j]) get multiplicity from dlist
351 * and store in multiplicity[j]
354 int j, Dih, i;
355 char name[4];
357 j = 0;
358 for (Dih = 0; (Dih < NONCHI+maxchi); Dih++)
360 for (i = 0; (i < nlist); i++)
362 std::strncpy(name, dlist[i].name, 3);
363 name[3] = '\0';
364 if (((Dih < edOmega) ) ||
365 ((Dih == edOmega) && (has_dihedral(edOmega, &(dlist[i])))) ||
366 ((Dih > edOmega) && (dlist[i].atm.Cn[Dih-NONCHI+3] != -1)))
368 /* default - we will correct the rest below */
369 multiplicity[j] = 3;
371 /* make omegas 2fold, though doesn't make much more sense than 3 */
372 if (Dih == edOmega && (has_dihedral(edOmega, &(dlist[i]))))
374 multiplicity[j] = 2;
377 /* dihedrals to aromatic rings, COO, CONH2 or guanidinium are 2fold*/
378 if (Dih > edOmega && (dlist[i].atm.Cn[Dih-NONCHI+3] != -1))
380 if ( ((std::strstr(name, "PHE") != nullptr) && (Dih == edChi2)) ||
381 ((std::strstr(name, "TYR") != nullptr) && (Dih == edChi2)) ||
382 ((std::strstr(name, "PTR") != nullptr) && (Dih == edChi2)) ||
383 ((std::strstr(name, "TRP") != nullptr) && (Dih == edChi2)) ||
384 ((std::strstr(name, "HIS") != nullptr) && (Dih == edChi2)) ||
385 ((std::strstr(name, "GLU") != nullptr) && (Dih == edChi3)) ||
386 ((std::strstr(name, "ASP") != nullptr) && (Dih == edChi2)) ||
387 ((std::strstr(name, "GLN") != nullptr) && (Dih == edChi3)) ||
388 ((std::strstr(name, "ASN") != nullptr) && (Dih == edChi2)) ||
389 ((std::strstr(name, "ARG") != nullptr) && (Dih == edChi4)) )
391 multiplicity[j] = 2;
394 j++;
398 if (j < nangles)
400 fprintf(stderr, "WARNING: not all dihedrals found in topology (only %d out of %d)!\n",
401 j, nangles);
403 /* Check for remaining dihedrals */
404 for (; (j < nangles); j++)
406 multiplicity[j] = 3;
411 void mk_chi_lookup (int **lookup, int maxchi,
412 int nlist, t_dlist dlist[])
415 /* by grs. should rewrite everything to use this. (but haven't,
416 * and at mmt only used in get_chi_product_traj
417 * returns the dihed number given the residue number (from-0)
418 * and chi (from-0) nr. -1 for chi undefined for that res (eg gly, ala..)*/
420 int i, j, Dih, Chi;
422 j = 0;
423 /* NONCHI points to chi1, therefore we have to start counting there. */
424 for (Dih = NONCHI; (Dih < NONCHI+maxchi); Dih++)
426 for (i = 0; (i < nlist); i++)
428 Chi = Dih - NONCHI;
429 if (((Dih < edOmega) ) ||
430 ((Dih == edOmega) && (has_dihedral(edOmega, &(dlist[i])))) ||
431 ((Dih > edOmega) && (dlist[i].atm.Cn[Dih-NONCHI+3] != -1)))
433 /* grs debug printf("Not OK? i %d j %d Dih %d \n", i, j, Dih) ; */
434 if (Dih > edOmega)
436 lookup[i][Chi] = j;
438 j++;
440 else
442 lookup[i][Chi] = -1;
450 void get_chi_product_traj (real **dih, int nframes, int nlist,
451 int maxchi, t_dlist dlist[], real time[],
452 int **lookup, int *multiplicity, gmx_bool bRb, gmx_bool bNormalize,
453 real core_frac, gmx_bool bAll, const char *fnall,
454 const gmx_output_env_t *oenv)
457 gmx_bool bRotZero, bHaveChi = FALSE;
458 int accum = 0, index, i, j, k, Xi, n, b;
459 real *chi_prtrj;
460 int *chi_prhist;
461 int nbin;
462 FILE *fp, *fpall;
463 char hisfile[256], histitle[256], *namept;
465 int (*calc_bin)(real, int, real);
467 /* Analysis of dihedral transitions */
468 fprintf(stderr, "Now calculating Chi product trajectories...\n");
470 if (bRb)
472 calc_bin = calc_RBbin;
474 else
476 calc_bin = calc_Nbin;
479 snew(chi_prtrj, nframes);
481 /* file for info on all residues */
482 if (bNormalize)
484 fpall = xvgropen(fnall, "Cumulative Rotamers", "Residue", "Probability", oenv);
486 else
488 fpall = xvgropen(fnall, "Cumulative Rotamers", "Residue", "# Counts", oenv);
491 for (i = 0; (i < nlist); i++)
494 /* get nbin, the nr. of cumulative rotamers that need to be considered */
495 nbin = 1;
496 for (Xi = 0; Xi < maxchi; Xi++)
498 index = lookup[i][Xi]; /* chi_(Xi+1) of res i (-1 if off end) */
499 if (index >= 0)
501 n = multiplicity[index];
502 nbin = n*nbin;
505 nbin += 1; /* for the "zero rotamer", outside the core region */
507 for (j = 0; (j < nframes); j++)
510 bRotZero = FALSE;
511 bHaveChi = TRUE;
512 index = lookup[i][0]; /* index into dih of chi1 of res i */
513 if (index == -1)
515 bRotZero = TRUE;
516 bHaveChi = FALSE;
518 else
520 b = calc_bin(dih[index][j], multiplicity[index], core_frac);
521 accum = b - 1;
522 if (b == 0)
524 bRotZero = TRUE;
526 for (Xi = 1; Xi < maxchi; Xi++)
528 index = lookup[i][Xi]; /* chi_(Xi+1) of res i (-1 if off end) */
529 if (index >= 0)
531 n = multiplicity[index];
532 b = calc_bin(dih[index][j], n, core_frac);
533 accum = n * accum + b - 1;
534 if (b == 0)
536 bRotZero = TRUE;
540 accum++;
542 if (bRotZero)
544 chi_prtrj[j] = 0.0;
546 else
548 chi_prtrj[j] = accum;
549 if (accum+1 > nbin)
551 nbin = accum+1;
555 if (bHaveChi)
558 if (bAll)
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(nullptr, nframes, chi_prtrj, nbin, chi_prhist, 0, nbin);
568 if (bAll)
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++)
582 if (bNormalize)
584 fprintf(fp, "%5d %10g\n", k, (1.0*chi_prhist[k])/nframes);
586 else
588 fprintf(fp, "%5d %10d\n", k, chi_prhist[k]);
591 fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
592 xvgrclose(fp);
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++)
602 if (bNormalize)
604 fprintf(fpall, " %10g", (1.0*chi_prhist[k])/nframes);
606 else
608 fprintf(fpall, " %10d", chi_prhist[k]);
611 fprintf(fpall, "\n");
613 sfree(chi_prhist);
614 /* histogram done */
618 sfree(chi_prtrj);
619 xvgrclose(fpall);
620 fprintf(stderr, "\n");
624 void calc_distribution_props(int nh, int histo[], real start,
625 int nkkk, t_karplus kkk[],
626 real *S2)
628 real d, dc, ds, c1, c2, tdc, tds;
629 real fac, ang, invth, Jc;
630 int i, j, th;
632 if (nh == 0)
634 gmx_fatal(FARGS, "No points in histogram (%s, %d)", __FILE__, __LINE__);
636 fac = 2*M_PI/nh;
638 /* Compute normalisation factor */
639 th = 0;
640 for (j = 0; (j < nh); j++)
642 th += histo[j];
644 invth = 1.0/th;
646 for (i = 0; (i < nkkk); i++)
648 kkk[i].Jc = 0;
649 kkk[i].Jcsig = 0;
651 tdc = 0, tds = 0;
652 for (j = 0; (j < nh); j++)
654 d = invth*histo[j];
655 ang = j*fac-start;
656 c1 = std::cos(ang);
657 dc = d*c1;
658 ds = d*std::sin(ang);
659 tdc += dc;
660 tds += ds;
661 for (i = 0; (i < nkkk); i++)
663 c1 = std::cos(ang+kkk[i].offset);
664 c2 = c1*c1;
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]*gmx::square(Jc);
670 for (i = 0; (i < nkkk); i++)
672 kkk[i].Jc /= th;
673 kkk[i].Jcsig = std::sqrt(kkk[i].Jcsig/th-gmx::square(kkk[i].Jc));
675 *S2 = tdc*tdc+tds*tds;
678 static void calc_angles(struct t_pbc *pbc,
679 int n3, int index[], real ang[], rvec x_s[])
681 int i, ix, t1, t2;
682 rvec r_ij, r_kj;
683 real costh = 0.0;
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);
690 if (debug)
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)
701 int i;
702 real trans = 0, gauche = 0;
703 real angle;
705 for (i = 0; i < nangles; i++)
707 angle = angles[i] * RAD2DEG;
709 if (angle > 135 && angle < 225)
711 trans += 1.0;
713 else if (angle > 270 && angle < 330)
715 gauche += 1.0;
717 else if (angle < 90 && angle > 30)
719 gauche += 1.0;
722 if (trans+gauche > 0)
724 return trans/(trans+gauche);
726 else
728 return 0;
732 static void calc_dihs(struct t_pbc *pbc,
733 int n4, int index[], real ang[], rvec x_s[])
735 int i, ix, t1, t2, t3;
736 rvec r_ij, r_kj, r_kl, m, n;
737 real sign, aaa;
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)
754 double dx;
755 int i, ind;
757 if (minx == maxx)
759 minx = maxx = data[0];
760 for (i = 1; (i < ndata); i++)
762 minx = std::min(minx, data[i]);
763 maxx = std::max(maxx, data[i]);
765 fprintf(log, "Min data: %10g Max data: %10g\n", minx, maxx);
767 dx = npoints/(maxx-minx);
768 if (debug)
770 fprintf(debug,
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 = static_cast<int>((data[i]-minx)*dx);
777 if ((ind >= 0) && (ind < npoints))
779 histo[ind]++;
781 else
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[])
790 int i;
791 double d, fac;
793 d = 0;
794 for (i = 0; (i < npoints); i++)
796 d += dx*histo[i];
798 if (d == 0)
800 fprintf(stderr, "Empty histogram!\n");
801 return;
803 fac = 1.0/d;
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, int index[],
815 real **trans_frac,
816 real **aver_angle,
817 real *dih[],
818 const gmx_output_env_t *oenv)
820 struct t_pbc *pbc;
821 t_trxstatus *status;
822 int i, angind, total, teller;
823 int nangles, n_alloc;
824 real t, fraction, pifac, aa, angle;
825 real *angles[2];
826 matrix box;
827 rvec *x;
828 int cur = 0;
829 #define prev (1-cur)
831 snew(pbc, 1);
832 read_first_x(oenv, &status, trj_fn, &t, &x, box);
834 if (bAngles)
836 nangles = isize/3;
837 pifac = M_PI;
839 else
841 nangles = isize/4;
842 pifac = 2.0*M_PI;
844 snew(angles[cur], nangles);
845 snew(angles[prev], nangles);
847 /* Start the loop over frames */
848 total = 0;
849 teller = 0;
850 n_alloc = 0;
851 *time = nullptr;
852 *trans_frac = nullptr;
853 *aver_angle = nullptr;
857 if (teller >= n_alloc)
859 n_alloc += 100;
860 if (bSaveAll)
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);
872 (*time)[teller] = t;
874 if (pbc)
876 set_pbc(pbc, -1, box);
879 if (bAngles)
881 calc_angles(pbc, isize, index, angles[cur], x);
883 else
885 calc_dihs(pbc, isize, index, angles[cur], x);
887 /* Trans fraction */
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].
897 if (bRb)
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... */
909 if (bPBC)
911 for (i = 0; (i < nangles); i++)
913 real dd = angles[cur][i];
914 angles[cur][i] = std::atan2(std::sin(dd), std::cos(dd));
917 else
919 if (teller > 1)
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;
936 /* Average angles */
937 aa = 0;
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
945 -Pi when plotting.
948 angle = angles[cur][i];
949 if (!bAngles)
951 while (angle < -M_PI)
953 angle += 2*M_PI;
955 while (angle >= M_PI)
957 angle -= 2*M_PI;
960 angle += M_PI;
963 /* Update the distribution histogram */
964 angind = static_cast<int>((angle*maxangstat)/pifac + 0.5);
965 if (angind == maxangstat)
967 angind = 0;
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);
976 angstat[angind]++;
977 if (angind == maxangstat)
979 fprintf(stderr, "angle %d fr %d = %g\n", i, cur, angle);
982 total++;
985 /* average over all angles */
986 (*aver_angle)[teller] = (aa/nangles);
988 /* this copies all current dih. angles to dih[i], teller is frame */
989 if (bSaveAll)
991 for (i = 0; i < nangles; i++)
993 dih[i][teller] = angles[cur][i];
997 /* Swap buffers */
998 cur = prev;
1000 /* Increment loop counter */
1001 teller++;
1003 while (read_next_x(oenv, status, &t, x, box));
1004 close_trx(status);
1006 sfree(x);
1007 sfree(angles[cur]);
1008 sfree(angles[prev]);
1010 *nframes = teller;