Avoid writing xvgr formatting with -xvg none
[alexxy/gromacs.git] / src / tools / anadih.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team,
6  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012,2013, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <math.h>
43 #include <stdio.h>
44 #include <string.h>
45 #include "physics.h"
46 #include "smalloc.h"
47 #include "macros.h"
48 #include "txtdump.h"
49 #include "bondf.h"
50 #include "xvgr.h"
51 #include "typedefs.h"
52 #include "vec.h"
53 #include "gstat.h"
54 #include "confio.h"
55
56 void print_one(const output_env_t oenv, const char *base, const char *name,
57                const char *title, const char *ylabel, int nf, real time[],
58                real data[])
59 {
60     FILE *fp;
61     char  buf[256], t2[256];
62     int   k;
63
64     sprintf(buf, "%s%s.xvg", base, name);
65     fprintf(stderr, "\rPrinting %s  ", buf);
66     sprintf(t2, "%s %s", title, name);
67     fp = xvgropen(buf, t2, "Time (ps)", ylabel, oenv);
68     for (k = 0; (k < nf); k++)
69     {
70         fprintf(fp, "%10g  %10g\n", time[k], data[k]);
71     }
72     ffclose(fp);
73 }
74
75 static int calc_RBbin(real phi, int multiplicity, real core_frac)
76 {
77     /* multiplicity and core_frac NOT used,
78      * just given to enable use of pt-to-fn in caller low_ana_dih_trans*/
79     static const real r30  = M_PI/6.0;
80     static const real r90  = M_PI/2.0;
81     static const real r150 = M_PI*5.0/6.0;
82
83     if ((phi < r30) && (phi > -r30))
84     {
85         return 1;
86     }
87     else if ((phi > -r150) && (phi < -r90))
88     {
89         return 2;
90     }
91     else if ((phi < r150) && (phi > r90))
92     {
93         return 3;
94     }
95     return 0;
96 }
97
98 static int calc_Nbin(real phi, int multiplicity, real core_frac)
99 {
100     static const real r360 = 360*DEG2RAD;
101     real              rot_width, core_width, core_offset, low, hi;
102     int               bin;
103     /* with multiplicity 3 and core_frac 0.5
104      * 0<g(-)<120, 120<t<240, 240<g(+)<360
105      * 0< bin0 < 30, 30<bin1<90, 90<bin0<150, 150<bin2<210, 210<bin0<270, 270<bin3<330, 330<bin0<360
106      * so with multiplicity 3, bin1 is core g(-), bin2 is core t, bin3 is
107        core g(+), bin0 is between rotamers */
108     if (phi < 0)
109     {
110         phi += r360;
111     }
112
113     rot_width   = 360/multiplicity;
114     core_width  = core_frac * rot_width;
115     core_offset = (rot_width - core_width)/2.0;
116     for (bin = 1; bin <= multiplicity; bin++)
117     {
118         low  = ((bin - 1) * rot_width ) + core_offset;
119         hi   = ((bin - 1) * rot_width ) + core_offset + core_width;
120         low *= DEG2RAD;
121         hi  *= DEG2RAD;
122         if ((phi > low) && (phi < hi))
123         {
124             return bin;
125         }
126     }
127     return 0;
128 }
129
130 void ana_dih_trans(const char *fn_trans, const char *fn_histo,
131                    real **dih, int nframes, int nangles,
132                    const char *grpname, real *time, gmx_bool bRb,
133                    const output_env_t oenv)
134 {
135     /* just a wrapper; declare extra args, then chuck away at end. */
136     int      maxchi = 0;
137     t_dlist *dlist;
138     int     *multiplicity;
139     int      nlist = nangles;
140     int      k;
141
142     snew(dlist, nlist);
143     snew(multiplicity, nangles);
144     for (k = 0; (k < nangles); k++)
145     {
146         multiplicity[k] = 3;
147     }
148
149     low_ana_dih_trans(TRUE, fn_trans, TRUE, fn_histo, maxchi,
150                       dih, nlist, dlist, nframes,
151                       nangles, grpname, multiplicity, time, bRb, 0.5, oenv);
152     sfree(dlist);
153     sfree(multiplicity);
154
155 }
156
157 void low_ana_dih_trans(gmx_bool bTrans, const char *fn_trans,
158                        gmx_bool bHisto, const char *fn_histo, int maxchi,
159                        real **dih, int nlist, t_dlist dlist[], int nframes,
160                        int nangles, const char *grpname, int multiplicity[],
161                        real *time, gmx_bool bRb, real core_frac,
162                        const output_env_t oenv)
163 {
164     FILE *fp;
165     int  *tr_f, *tr_h;
166     char  title[256];
167     int   i, j, k, Dih, ntrans;
168     int   cur_bin, new_bin;
169     real  ttime, tt;
170     real *rot_occ[NROT];
171     int   (*calc_bin)(real, int, real);
172     real  dt;
173
174     if (1 <= nframes)
175     {
176         return;
177     }
178     /* Assumes the frames are equally spaced in time */
179     dt = (time[nframes-1]-time[0])/(nframes-1);
180
181     /* Analysis of dihedral transitions */
182     fprintf(stderr, "Now calculating transitions...\n");
183
184     if (bRb)
185     {
186         calc_bin = calc_RBbin;
187     }
188     else
189     {
190         calc_bin = calc_Nbin;
191     }
192
193     for (k = 0; k < NROT; k++)
194     {
195         snew(rot_occ[k], nangles);
196         for (i = 0; (i < nangles); i++)
197         {
198             rot_occ[k][i] = 0;
199         }
200     }
201     snew(tr_h, nangles);
202     snew(tr_f, nframes);
203
204     /* dih[i][j] is the dihedral angle i in frame j  */
205     ntrans = 0;
206     for (i = 0; (i < nangles); i++)
207     {
208
209         /*#define OLDIE*/
210 #ifdef OLDIE
211         mind = maxd = prev = dih[i][0];
212 #else
213         cur_bin = calc_bin(dih[i][0], multiplicity[i], core_frac);
214         rot_occ[cur_bin][i]++;
215 #endif
216         for (j = 1; (j < nframes); j++)
217         {
218             new_bin = calc_bin(dih[i][j], multiplicity[i], core_frac);
219             rot_occ[new_bin][i]++;
220 #ifndef OLDIE
221             if (cur_bin == 0)
222             {
223                 cur_bin = new_bin;
224             }
225             else if ((new_bin != 0) && (cur_bin != new_bin))
226             {
227                 cur_bin = new_bin;
228                 tr_f[j]++;
229                 tr_h[i]++;
230                 ntrans++;
231             }
232 #else
233             /* why is all this md rubbish periodic? Remove 360 degree periodicity */
234             if ( (dih[i][j] - prev) > M_PI)
235             {
236                 dih[i][j] -= 2*M_PI;
237             }
238             else if ( (dih[i][j] - prev) < -M_PI)
239             {
240                 dih[i][j] += 2*M_PI;
241             }
242
243             prev = dih[i][j];
244
245             mind = min(mind, dih[i][j]);
246             maxd = max(maxd, dih[i][j]);
247             if ( (maxd - mind) > 2*M_PI/3) /* or 120 degrees, assuming       */
248             {                              /* multiplicity 3. Not so general.*/
249                 tr_f[j]++;
250                 tr_h[i]++;
251                 maxd = mind = dih[i][j]; /* get ready for next transition  */
252                 ntrans++;
253             }
254 #endif
255         } /* end j */
256         for (k = 0; k < NROT; k++)
257         {
258             rot_occ[k][i] /= nframes;
259         }
260     } /* end i */
261     fprintf(stderr, "Total number of transitions: %10d\n", ntrans);
262     if (ntrans > 0)
263     {
264         ttime = (dt*nframes*nangles)/ntrans;
265         fprintf(stderr, "Time between transitions:    %10.3f ps\n", ttime);
266     }
267
268     /* new by grs - copy transitions from tr_h[] to dlist->ntr[]
269      * and rotamer populations from rot_occ to dlist->rot_occ[]
270      * based on fn histogramming in g_chi. diff roles for i and j here */
271
272     j = 0;
273     for (Dih = 0; (Dih < NONCHI+maxchi); Dih++)
274     {
275         for (i = 0; (i < nlist); i++)
276         {
277             if (((Dih  < edOmega) ) ||
278                 ((Dih == edOmega) && (has_dihedral(edOmega, &(dlist[i])))) ||
279                 ((Dih  > edOmega) && (dlist[i].atm.Cn[Dih-NONCHI+3] != -1)))
280             {
281                 /* grs debug  printf("Not OK? i %d j %d Dih %d \n", i, j, Dih) ; */
282                 dlist[i].ntr[Dih] = tr_h[j];
283                 for (k = 0; k < NROT; k++)
284                 {
285                     dlist[i].rot_occ[Dih][k] = rot_occ[k][j];
286                 }
287                 j++;
288             }
289         }
290     }
291
292     /* end addition by grs */
293
294     if (bTrans)
295     {
296         sprintf(title, "Number of transitions: %s", grpname);
297         fp = xvgropen(fn_trans, title, "Time (ps)", "# transitions/timeframe", oenv);
298         for (j = 0; (j < nframes); j++)
299         {
300             fprintf(fp, "%10.3f  %10d\n", time[j], tr_f[j]);
301         }
302         ffclose(fp);
303     }
304
305     /* Compute histogram from # transitions per dihedral */
306     /* Use old array */
307     for (j = 0; (j < nframes); j++)
308     {
309         tr_f[j] = 0;
310     }
311     for (i = 0; (i < nangles); i++)
312     {
313         tr_f[tr_h[i]]++;
314     }
315     for (j = nframes; ((tr_f[j-1] == 0) && (j > 0)); j--)
316     {
317         ;
318     }
319
320     ttime = dt*nframes;
321     if (bHisto)
322     {
323         sprintf(title, "Transition time: %s", grpname);
324         fp = xvgropen(fn_histo, title, "Time (ps)", "#", oenv);
325         for (i = j-1; (i > 0); i--)
326         {
327             if (tr_f[i] != 0)
328             {
329                 fprintf(fp, "%10.3f  %10d\n", ttime/i, tr_f[i]);
330             }
331         }
332         ffclose(fp);
333     }
334
335     sfree(tr_f);
336     sfree(tr_h);
337     for (k = 0; k < NROT; k++)
338     {
339         sfree(rot_occ[k]);
340     }
341
342 }
343
344 void mk_multiplicity_lookup (int *multiplicity, int maxchi, real **dih,
345                              int nlist, t_dlist dlist[], int nangles)
346 {
347     /* new by grs - for dihedral j (as in dih[j]) get multiplicity from dlist
348      * and store in multiplicity[j]
349      */
350
351     int  j, Dih, i;
352     char name[4];
353
354     j = 0;
355     for (Dih = 0; (Dih < NONCHI+maxchi); Dih++)
356     {
357         for (i = 0; (i < nlist); i++)
358         {
359             strncpy(name, dlist[i].name, 3);
360             name[3] = '\0';
361             if (((Dih  < edOmega) ) ||
362                 ((Dih == edOmega) && (has_dihedral(edOmega, &(dlist[i])))) ||
363                 ((Dih  > edOmega) && (dlist[i].atm.Cn[Dih-NONCHI+3] != -1)))
364             {
365                 /* default - we will correct the rest below */
366                 multiplicity[j] = 3;
367
368                 /* make omegas 2fold, though doesn't make much more sense than 3 */
369                 if (Dih == edOmega && (has_dihedral(edOmega, &(dlist[i]))))
370                 {
371                     multiplicity[j] = 2;
372                 }
373
374                 /* dihedrals to aromatic rings, COO, CONH2 or guanidinium are 2fold*/
375                 if (Dih > edOmega && (dlist[i].atm.Cn[Dih-NONCHI+3] != -1))
376                 {
377                     if ( ((strstr(name, "PHE") != NULL) && (Dih == edChi2))  ||
378                          ((strstr(name, "TYR") != NULL) && (Dih == edChi2))  ||
379                          ((strstr(name, "PTR") != NULL) && (Dih == edChi2))  ||
380                          ((strstr(name, "TRP") != NULL) && (Dih == edChi2))  ||
381                          ((strstr(name, "HIS") != NULL) && (Dih == edChi2))  ||
382                          ((strstr(name, "GLU") != NULL) && (Dih == edChi3))  ||
383                          ((strstr(name, "ASP") != NULL) && (Dih == edChi2))  ||
384                          ((strstr(name, "GLN") != NULL) && (Dih == edChi3))  ||
385                          ((strstr(name, "ASN") != NULL) && (Dih == edChi2))  ||
386                          ((strstr(name, "ARG") != NULL) && (Dih == edChi4))  )
387                     {
388                         multiplicity[j] = 2;
389                     }
390                 }
391                 j++;
392             }
393         }
394     }
395     if (j < nangles)
396     {
397         fprintf(stderr, "WARNING: not all dihedrals found in topology (only %d out of %d)!\n",
398                 j, nangles);
399     }
400     /* Check for remaining dihedrals */
401     for (; (j < nangles); j++)
402     {
403         multiplicity[j] = 3;
404     }
405
406 }
407
408 void mk_chi_lookup (int **lookup, int maxchi, real **dih,
409                     int nlist, t_dlist dlist[])
410 {
411
412     /* by grs. should rewrite everything to use this. (but haven't,
413      * and at mmt only used in get_chi_product_traj
414      * returns the dihed number given the residue number (from-0)
415      * and chi (from-0) nr. -1 for chi undefined for that res (eg gly, ala..)*/
416
417     int i, j, Dih, Chi;
418
419     j = 0;
420     for (Dih = 0; (Dih < NONCHI+maxchi); Dih++)
421     {
422         for (i = 0; (i < nlist); i++)
423         {
424             Chi = Dih - NONCHI;
425             if (((Dih  < edOmega) ) ||
426                 ((Dih == edOmega) && (has_dihedral(edOmega, &(dlist[i])))) ||
427                 ((Dih  > edOmega) && (dlist[i].atm.Cn[Dih-NONCHI+3] != -1)))
428             {
429                 /* grs debug  printf("Not OK? i %d j %d Dih %d \n", i, j, Dih) ; */
430                 if (Dih > edOmega)
431                 {
432                     lookup[i][Chi] = j;
433                 }
434                 j++;
435             }
436             else
437             {
438                 lookup[i][Chi] = -1;
439             }
440         }
441     }
442
443 }
444
445
446 void get_chi_product_traj (real **dih, int nframes, int nangles, int nlist,
447                            int maxchi, t_dlist dlist[], real time[],
448                            int **lookup, int *multiplicity, gmx_bool bRb, gmx_bool bNormalize,
449                            real core_frac, gmx_bool bAll, const char *fnall,
450                            const output_env_t oenv)
451 {
452
453     gmx_bool bRotZero, bHaveChi = FALSE;
454     int      accum = 0, index, i, j, k, Xi, n, b;
455     real    *chi_prtrj;
456     int     *chi_prhist;
457     int      nbin;
458     FILE    *fp, *fpall;
459     char     hisfile[256], histitle[256], *namept;
460
461     int      (*calc_bin)(real, int, real);
462
463     /* Analysis of dihedral transitions */
464     fprintf(stderr, "Now calculating Chi product trajectories...\n");
465
466     if (bRb)
467     {
468         calc_bin = calc_RBbin;
469     }
470     else
471     {
472         calc_bin = calc_Nbin;
473     }
474
475     snew(chi_prtrj, nframes);
476
477     /* file for info on all residues */
478     if (bNormalize)
479     {
480         fpall = xvgropen(fnall, "Cumulative Rotamers", "Residue", "Probability", oenv);
481     }
482     else
483     {
484         fpall = xvgropen(fnall, "Cumulative Rotamers", "Residue", "# Counts", oenv);
485     }
486
487     for (i = 0; (i < nlist); i++)
488     {
489
490         /* get nbin, the nr. of cumulative rotamers that need to be considered */
491         nbin = 1;
492         for (Xi = 0; Xi < maxchi; Xi++)
493         {
494             index = lookup[i][Xi]; /* chi_(Xi+1) of res i (-1 if off end) */
495             if (index >= 0)
496             {
497                 n    = multiplicity[index];
498                 nbin = n*nbin;
499             }
500         }
501         nbin += 1; /* for the "zero rotamer", outside the core region */
502
503         for (j = 0; (j < nframes); j++)
504         {
505
506             bRotZero = FALSE;
507             bHaveChi = TRUE;
508             index    = lookup[i][0]; /* index into dih of chi1 of res i */
509             if (index == -1)
510             {
511                 b        = 0;
512                 bRotZero = TRUE;
513                 bHaveChi = FALSE;
514             }
515             else
516             {
517                 b     = calc_bin(dih[index][j], multiplicity[index], core_frac);
518                 accum = b - 1;
519                 if (b == 0)
520                 {
521                     bRotZero = TRUE;
522                 }
523                 for (Xi = 1; Xi < maxchi; Xi++)
524                 {
525                     index = lookup[i][Xi]; /* chi_(Xi+1) of res i (-1 if off end) */
526                     if (index >= 0)
527                     {
528                         n     = multiplicity[index];
529                         b     = calc_bin(dih[index][j], n, core_frac);
530                         accum = n * accum + b - 1;
531                         if (b == 0)
532                         {
533                             bRotZero = TRUE;
534                         }
535                     }
536                 }
537                 accum++;
538             }
539             if (bRotZero)
540             {
541                 chi_prtrj[j] = 0.0;
542             }
543             else
544             {
545                 chi_prtrj[j] = accum;
546                 if (accum+1 > nbin)
547                 {
548                     nbin = accum+1;
549                 }
550             }
551         }
552         if (bHaveChi)
553         {
554
555             if (bAll)
556             {
557                 /* print cuml rotamer vs time */
558                 print_one(oenv, "chiproduct", dlist[i].name, "chi product for",
559                           "cumulative rotamer", nframes, time, chi_prtrj);
560             }
561
562             /* make a histogram pf culm. rotamer occupancy too */
563             snew(chi_prhist, nbin);
564             make_histo(NULL, nframes, chi_prtrj, nbin, chi_prhist, 0, nbin);
565             if (bAll)
566             {
567                 sprintf(hisfile, "histo-chiprod%s.xvg", dlist[i].name);
568                 sprintf(histitle, "cumulative rotamer distribution for %s", dlist[i].name);
569                 fprintf(stderr, "  and %s  ", hisfile);
570                 fp = xvgropen(hisfile, histitle, "number", "", oenv);
571                 if(output_env_get_print_xvgr_codes(oenv))
572                 {
573                     fprintf(fp, "@ xaxis tick on\n");
574                     fprintf(fp, "@ xaxis tick major 1\n");
575                     fprintf(fp, "@ type xy\n");
576                 }
577                 for (k = 0; (k < nbin); k++)
578                 {
579                     if (bNormalize)
580                     {
581                         fprintf(fp, "%5d  %10g\n", k, (1.0*chi_prhist[k])/nframes);
582                     }
583                     else
584                     {
585                         fprintf(fp, "%5d  %10d\n", k, chi_prhist[k]);
586                     }
587                 }
588                 fprintf(fp, "%s\n", output_env_get_print_xvgr_codes(oenv) ? "&" : "");
589                 ffclose(fp);
590             }
591
592             /* and finally print out occupancies to a single file */
593             /* get the gmx from-1 res nr by setting a ptr to the number part
594              * of dlist[i].name - potential bug for 4-letter res names... */
595             namept = dlist[i].name + 3;
596             fprintf(fpall, "%5s ", namept);
597             for (k = 0; (k < nbin); k++)
598             {
599                 if (bNormalize)
600                 {
601                     fprintf(fpall, "  %10g", (1.0*chi_prhist[k])/nframes);
602                 }
603                 else
604                 {
605                     fprintf(fpall, "  %10d", chi_prhist[k]);
606                 }
607             }
608             fprintf(fpall, "\n");
609
610             sfree(chi_prhist);
611             /* histogram done */
612         }
613     }
614
615     sfree(chi_prtrj);
616     ffclose(fpall);
617     fprintf(stderr, "\n");
618
619 }
620
621 void calc_distribution_props(int nh, int histo[], real start,
622                              int nkkk, t_karplus kkk[],
623                              real *S2)
624 {
625     real d, dc, ds, c1, c2, tdc, tds;
626     real fac, ang, invth, Jc;
627     int  i, j, th;
628
629     if (nh == 0)
630     {
631         gmx_fatal(FARGS, "No points in histogram (%s, %d)", __FILE__, __LINE__);
632     }
633     fac = 2*M_PI/nh;
634
635     /* Compute normalisation factor */
636     th = 0;
637     for (j = 0; (j < nh); j++)
638     {
639         th += histo[j];
640     }
641     invth = 1.0/th;
642
643     for (i = 0; (i < nkkk); i++)
644     {
645         kkk[i].Jc    = 0;
646         kkk[i].Jcsig = 0;
647     }
648     tdc = 0, tds = 0;
649     for (j = 0; (j < nh); j++)
650     {
651         d    = invth*histo[j];
652         ang  = j*fac-start;
653         c1   = cos(ang);
654         c2   = c1*c1;
655         dc   = d*c1;
656         ds   = d*sin(ang);
657         tdc += dc;
658         tds += ds;
659         for (i = 0; (i < nkkk); i++)
660         {
661             c1            = cos(ang+kkk[i].offset);
662             c2            = c1*c1;
663             Jc            = (kkk[i].A*c2 + kkk[i].B*c1 + kkk[i].C);
664             kkk[i].Jc    += histo[j]*Jc;
665             kkk[i].Jcsig += histo[j]*sqr(Jc);
666         }
667     }
668     for (i = 0; (i < nkkk); i++)
669     {
670         kkk[i].Jc    /= th;
671         kkk[i].Jcsig  = sqrt(kkk[i].Jcsig/th-sqr(kkk[i].Jc));
672     }
673     *S2 = tdc*tdc+tds*tds;
674 }
675
676 static void calc_angles(FILE *log, t_pbc *pbc,
677                         int n3, atom_id index[], real ang[], rvec x_s[])
678 {
679     int  i, ix, t1, t2;
680     rvec r_ij, r_kj;
681     real costh;
682
683     for (i = ix = 0; (ix < n3); i++, ix += 3)
684     {
685         ang[i] = bond_angle(x_s[index[ix]], x_s[index[ix+1]], x_s[index[ix+2]],
686                             pbc, r_ij, r_kj, &costh, &t1, &t2);
687     }
688     if (debug)
689     {
690         fprintf(debug, "Angle[0]=%g, costh=%g, index0 = %d, %d, %d\n",
691                 ang[0], costh, index[0], index[1], index[2]);
692         pr_rvec(debug, 0, "rij", r_ij, DIM, TRUE);
693         pr_rvec(debug, 0, "rkj", r_kj, DIM, TRUE);
694     }
695 }
696
697 static real calc_fraction(real angles[], int nangles)
698 {
699     int  i;
700     real trans = 0, gauche = 0;
701     real angle;
702
703     for (i = 0; i < nangles; i++)
704     {
705         angle = angles[i] * RAD2DEG;
706
707         if (angle > 135 && angle < 225)
708         {
709             trans += 1.0;
710         }
711         else if (angle > 270 && angle < 330)
712         {
713             gauche += 1.0;
714         }
715         else if (angle < 90 && angle > 30)
716         {
717             gauche += 1.0;
718         }
719     }
720     if (trans+gauche > 0)
721     {
722         return trans/(trans+gauche);
723     }
724     else
725     {
726         return 0;
727     }
728 }
729
730 static void calc_dihs(FILE *log, t_pbc *pbc,
731                       int n4, atom_id index[], real ang[], rvec x_s[])
732 {
733     int  i, ix, t1, t2, t3;
734     rvec r_ij, r_kj, r_kl, m, n;
735     real sign, aaa;
736
737     for (i = ix = 0; (ix < n4); i++, ix += 4)
738     {
739         aaa = dih_angle(x_s[index[ix]], x_s[index[ix+1]], x_s[index[ix+2]],
740                         x_s[index[ix+3]], pbc,
741                         r_ij, r_kj, r_kl, m, n,
742                         &sign, &t1, &t2, &t3);
743
744         ang[i] = aaa; /* not taking into account ryckaert bellemans yet */
745     }
746 }
747
748 void make_histo(FILE *log,
749                 int ndata, real data[], int npoints, int histo[],
750                 real minx, real maxx)
751 {
752     double dx;
753     int    i, ind;
754
755     if (minx == maxx)
756     {
757         minx = maxx = data[0];
758         for (i = 1; (i < ndata); i++)
759         {
760             minx = min(minx, data[i]);
761             maxx = max(maxx, data[i]);
762         }
763         fprintf(log, "Min data: %10g  Max data: %10g\n", minx, maxx);
764     }
765     dx = (double)npoints/(maxx-minx);
766     if (debug)
767     {
768         fprintf(debug,
769                 "Histogramming: ndata=%d, nhisto=%d, minx=%g,maxx=%g,dx=%g\n",
770                 ndata, npoints, minx, maxx, dx);
771     }
772     for (i = 0; (i < ndata); i++)
773     {
774         ind = (data[i]-minx)*dx;
775         if ((ind >= 0) && (ind < npoints))
776         {
777             histo[ind]++;
778         }
779         else
780         {
781             fprintf(log, "index = %d, data[%d] = %g\n", ind, i, data[i]);
782         }
783     }
784 }
785
786 void normalize_histo(int npoints, int histo[], real dx, real normhisto[])
787 {
788     int    i;
789     double d, fac;
790
791     d = 0;
792     for (i = 0; (i < npoints); i++)
793     {
794         d += dx*histo[i];
795     }
796     if (d == 0)
797     {
798         fprintf(stderr, "Empty histogram!\n");
799         return;
800     }
801     fac = 1.0/d;
802     for (i = 0; (i < npoints); i++)
803     {
804         normhisto[i] = fac*histo[i];
805     }
806 }
807
808 void read_ang_dih(const char *trj_fn,
809                   gmx_bool bAngles, gmx_bool bSaveAll, gmx_bool bRb, gmx_bool bPBC,
810                   int maxangstat, int angstat[],
811                   int *nframes, real **time,
812                   int isize, atom_id index[],
813                   real **trans_frac,
814                   real **aver_angle,
815                   real *dih[],
816                   const output_env_t oenv)
817 {
818     t_pbc       *pbc;
819     t_trxstatus *status;
820     int          i, angind, natoms, total, teller;
821     int          nangles, n_alloc;
822     real         t, fraction, pifac, aa, angle;
823     real        *angles[2];
824     matrix       box;
825     rvec        *x;
826     int          cur = 0;
827 #define prev (1-cur)
828
829     snew(pbc, 1);
830     natoms = read_first_x(oenv, &status, trj_fn, &t, &x, box);
831
832     if (bAngles)
833     {
834         nangles = isize/3;
835         pifac   = M_PI;
836     }
837     else
838     {
839         nangles = isize/4;
840         pifac   = 2.0*M_PI;
841     }
842     snew(angles[cur], nangles);
843     snew(angles[prev], nangles);
844
845     /* Start the loop over frames */
846     total       = 0;
847     teller      = 0;
848     n_alloc     = 0;
849     *time       = NULL;
850     *trans_frac = NULL;
851     *aver_angle = NULL;
852
853     do
854     {
855         if (teller >= n_alloc)
856         {
857             n_alloc += 100;
858             if (bSaveAll)
859             {
860                 for (i = 0; (i < nangles); i++)
861                 {
862                     srenew(dih[i], n_alloc);
863                 }
864             }
865             srenew(*time, n_alloc);
866             srenew(*trans_frac, n_alloc);
867             srenew(*aver_angle, n_alloc);
868         }
869
870         (*time)[teller] = t;
871
872         if (pbc)
873         {
874             set_pbc(pbc, -1, box);
875         }
876
877         if (bAngles)
878         {
879             calc_angles(stdout, pbc, isize, index, angles[cur], x);
880         }
881         else
882         {
883             calc_dihs(stdout, pbc, isize, index, angles[cur], x);
884
885             /* Trans fraction */
886             fraction              = calc_fraction(angles[cur], nangles);
887             (*trans_frac)[teller] = fraction;
888
889             /* Change Ryckaert-Bellemans dihedrals to polymer convention
890              * Modified 990913 by Erik:
891              * We actually shouldn't change the convention, since it's
892              * calculated from polymer above, but we change the intervall
893              * from [-180,180] to [0,360].
894              */
895             if (bRb)
896             {
897                 for (i = 0; (i < nangles); i++)
898                 {
899                     if (angles[cur][i] <= 0.0)
900                     {
901                         angles[cur][i] += 2*M_PI;
902                     }
903                 }
904             }
905
906             /* Periodicity in dihedral space... */
907             if (bPBC)
908             {
909                 for (i = 0; (i < nangles); i++)
910                 {
911                     real dd = angles[cur][i];
912                     angles[cur][i] = atan2(sin(dd), cos(dd));
913                 }
914             }
915             else
916             {
917                 if (teller > 1)
918                 {
919                     for (i = 0; (i < nangles); i++)
920                     {
921                         while (angles[cur][i] <= angles[prev][i] - M_PI)
922                         {
923                             angles[cur][i] += 2*M_PI;
924                         }
925                         while (angles[cur][i] > angles[prev][i] + M_PI)
926                         {
927                             angles[cur][i] -= 2*M_PI;
928                         }
929                     }
930                 }
931             }
932         }
933
934         /* Average angles */
935         aa = 0;
936         for (i = 0; (i < nangles); i++)
937         {
938             aa = aa+angles[cur][i];
939
940             /* angle in rad / 2Pi * max determines bin. bins go from 0 to maxangstat,
941                even though scale goes from -pi to pi (dihedral) or -pi/2 to pi/2
942                (angle) Basically: translate the x-axis by Pi. Translate it back by
943                -Pi when plotting.
944              */
945
946             angle = angles[cur][i];
947             if (!bAngles)
948             {
949                 while (angle < -M_PI)
950                 {
951                     angle += 2*M_PI;
952                 }
953                 while (angle >= M_PI)
954                 {
955                     angle -= 2*M_PI;
956                 }
957
958                 angle += M_PI;
959             }
960
961             /* Update the distribution histogram */
962             angind = (int) ((angle*maxangstat)/pifac + 0.5);
963             if (angind == maxangstat)
964             {
965                 angind = 0;
966             }
967             if ( (angind < 0) || (angind >= maxangstat) )
968             {
969                 /* this will never happen */
970                 gmx_fatal(FARGS, "angle (%f) index out of range (0..%d) : %d\n",
971                           angle, maxangstat, angind);
972             }
973
974             angstat[angind]++;
975             if (angind == maxangstat)
976             {
977                 fprintf(stderr, "angle %d fr %d = %g\n", i, cur, angle);
978             }
979
980             total++;
981         }
982
983         /* average over all angles */
984         (*aver_angle)[teller] = (aa/nangles);
985
986         /* this copies all current dih. angles to dih[i], teller is frame */
987         if (bSaveAll)
988         {
989             for (i = 0; i < nangles; i++)
990             {
991                 dih[i][teller] = angles[cur][i];
992             }
993         }
994
995         /* Swap buffers */
996         cur = prev;
997
998         /* Increment loop counter */
999         teller++;
1000     }
1001     while (read_next_x(oenv, status, &t, natoms, x, box));
1002     close_trj(status);
1003
1004     sfree(x);
1005     sfree(angles[cur]);
1006     sfree(angles[prev]);
1007
1008     *nframes = teller;
1009 }