int gmx_g_angle(int argc, char *argv[]) { static const char *desc[] = { "[TT]g_angle[tt] computes the angle distribution for a number of angles", "or dihedrals.[PAR]", "With option [TT]-ov[tt], you can plot the average angle of", "a group of angles as a function of time. With the [TT]-all[tt] option,", "the first graph is the average and the rest are the individual angles.[PAR]", "With the [TT]-of[tt] option, [TT]g_angle[tt] also calculates the fraction of trans", "dihedrals (only for dihedrals) as function of time, but this is", "probably only fun for a select few.[PAR]", "With option [TT]-oc[tt], a dihedral correlation function is calculated.[PAR]", "It should be noted that the index file must contain", "atom triplets for angles or atom quadruplets for dihedrals.", "If this is not the case, the program will crash.[PAR]", "With option [TT]-or[tt], a trajectory file is dumped containing cos and", "sin of selected dihedral angles, which subsequently can be used as", "input for a principal components analysis using [TT]g_covar[tt].[PAR]", "Option [TT]-ot[tt] plots when transitions occur between", "dihedral rotamers of multiplicity 3 and [TT]-oh[tt]", "records a histogram of the times between such transitions,", "assuming the input trajectory frames are equally spaced in time." }; static const char *opt[] = { NULL, "angle", "dihedral", "improper", "ryckaert-bellemans", NULL }; static gmx_bool bALL = FALSE, bChandler = FALSE, bAverCorr = FALSE, bPBC = TRUE; static real binwidth = 1; t_pargs pa[] = { { "-type", FALSE, etENUM, {opt}, "Type of angle to analyse" }, { "-all", FALSE, etBOOL, {&bALL}, "Plot all angles separately in the averages file, in the order of appearance in the index file." }, { "-binwidth", FALSE, etREAL, {&binwidth}, "binwidth (degrees) for calculating the distribution" }, { "-periodic", FALSE, etBOOL, {&bPBC}, "Print dihedral angles modulo 360 degrees" }, { "-chandler", FALSE, etBOOL, {&bChandler}, "Use Chandler correlation function (N[trans] = 1, N[gauche] = 0) rather than cosine correlation function. Trans is defined as phi < -60 or phi > 60." }, { "-avercorr", FALSE, etBOOL, {&bAverCorr}, "Average the correlation functions for the individual angles/dihedrals" } }; static const char *bugs[] = { "Counting transitions only works for dihedrals with multiplicity 3" }; FILE *out; real tmp, dt; int status, isize; atom_id *index; char *grpname; real maxang, Jc, S2, norm_fac, maxstat; unsigned long mode; int nframes, maxangstat, mult, *angstat; int i, j, total, nangles, natoms, nat2, first, last, angind; gmx_bool bAver, bRb, bPeriodic, bFrac, /* calculate fraction too? */ bTrans, /* worry about transtions too? */ bCorr; /* correlation function ? */ real t, aa, aver, aver2, aversig, fraction; /* fraction trans dihedrals */ double tfrac = 0; char title[256]; real **dih = NULL; /* mega array with all dih. angles at all times*/ char buf[80]; real *time, *trans_frac, *aver_angle; t_filenm fnm[] = { { efTRX, "-f", NULL, ffREAD }, { efNDX, NULL, "angle", ffREAD }, { efXVG, "-od", "angdist", ffWRITE }, { efXVG, "-ov", "angaver", ffOPTWR }, { efXVG, "-of", "dihfrac", ffOPTWR }, { efXVG, "-ot", "dihtrans", ffOPTWR }, { efXVG, "-oh", "trhisto", ffOPTWR }, { efXVG, "-oc", "dihcorr", ffOPTWR }, { efTRR, "-or", NULL, ffOPTWR } }; #define NFILE asize(fnm) int npargs; t_pargs *ppa; output_env_t oenv; npargs = asize(pa); ppa = add_acf_pargs(&npargs, pa); parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE, NFILE, fnm, npargs, ppa, asize(desc), desc, asize(bugs), bugs, &oenv); mult = 4; maxang = 360.0; bRb = FALSE; switch (opt[0][0]) { case 'a': mult = 3; maxang = 180.0; break; case 'd': break; case 'i': break; case 'r': bRb = TRUE; break; } if (opt2bSet("-or", NFILE, fnm)) { if (mult != 4) { gmx_fatal(FARGS, "Can not combine angles with trn dump"); } else { please_cite(stdout, "Mu2005a"); } } /* Calculate bin size */ maxangstat = (int)(maxang/binwidth+0.5); binwidth = maxang/maxangstat; rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &isize, &index, &grpname); nangles = isize/mult; if ((isize % mult) != 0) { gmx_fatal(FARGS, "number of index elements not multiple of %d, " "these can not be %s\n", mult, (mult == 3) ? "angle triplets" : "dihedral quadruplets"); } /* Check whether specific analysis has to be performed */ bCorr = opt2bSet("-oc", NFILE, fnm); bAver = opt2bSet("-ov", NFILE, fnm); bTrans = opt2bSet("-ot", NFILE, fnm); bFrac = opt2bSet("-of", NFILE, fnm); if (bTrans && opt[0][0] != 'd') { fprintf(stderr, "Option -ot should only accompany -type dihedral. Disabling -ot.\n"); bTrans = FALSE; } if (bChandler && !bCorr) { bCorr = TRUE; } if (bFrac && !bRb) { fprintf(stderr, "Warning:" " calculating fractions as defined in this program\n" "makes sense for Ryckaert Bellemans dihs. only. Ignoring -of\n\n"); bFrac = FALSE; } if ( (bTrans || bFrac || bCorr) && mult == 3) { gmx_fatal(FARGS, "Can only do transition, fraction or correlation\n" "on dihedrals. Select -d\n"); } /* * We need to know the nr of frames so we can allocate memory for an array * with all dihedral angles at all timesteps. Works for me. */ if (bTrans || bCorr || bALL || opt2bSet("-or", NFILE, fnm)) { snew(dih, nangles); } snew(angstat, maxangstat); read_ang_dih(ftp2fn(efTRX, NFILE, fnm), (mult == 3), bALL || bCorr || bTrans || opt2bSet("-or", NFILE, fnm), bRb, bPBC, maxangstat, angstat, &nframes, &time, isize, index, &trans_frac, &aver_angle, dih, oenv); dt = (time[nframes-1]-time[0])/(nframes-1); if (bAver) { sprintf(title, "Average Angle: %s", grpname); out = xvgropen(opt2fn("-ov", NFILE, fnm), title, "Time (ps)", "Angle (degrees)", oenv); for (i = 0; (i < nframes); i++) { fprintf(out, "%10.5f %8.3f", time[i], aver_angle[i]*RAD2DEG); if (bALL) { for (j = 0; (j < nangles); j++) { if (bPBC) { real dd = dih[j][i]; fprintf(out, " %8.3f", atan2(sin(dd), cos(dd))*RAD2DEG); } else { fprintf(out, " %8.3f", dih[j][i]*RAD2DEG); } } } fprintf(out, "\n"); } ffclose(out); } if (opt2bSet("-or", NFILE, fnm)) { dump_dih_trn(nframes, nangles, dih, opt2fn("-or", NFILE, fnm), time); } if (bFrac) { sprintf(title, "Trans fraction: %s", grpname); out = xvgropen(opt2fn("-of", NFILE, fnm), title, "Time (ps)", "Fraction", oenv); tfrac = 0.0; for (i = 0; (i < nframes); i++) { fprintf(out, "%10.5f %10.3f\n", time[i], trans_frac[i]); tfrac += trans_frac[i]; } ffclose(out); tfrac /= nframes; fprintf(stderr, "Average trans fraction: %g\n", tfrac); } sfree(trans_frac); if (bTrans) { ana_dih_trans(opt2fn("-ot", NFILE, fnm), opt2fn("-oh", NFILE, fnm), dih, nframes, nangles, grpname, time, bRb, oenv); } if (bCorr) { /* Autocorrelation function */ if (nframes < 2) { fprintf(stderr, "Not enough frames for correlation function\n"); } else { if (bChandler) { real dval, sixty = DEG2RAD*60; gmx_bool bTest; for (i = 0; (i < nangles); i++) { for (j = 0; (j < nframes); j++) { dval = dih[i][j]; if (bRb) { bTest = (dval > -sixty) && (dval < sixty); } else { bTest = (dval < -sixty) || (dval > sixty); } if (bTest) { dih[i][j] = dval-tfrac; } else { dih[i][j] = -tfrac; } } } } if (bChandler) { mode = eacNormal; } else { mode = eacCos; } do_autocorr(opt2fn("-oc", NFILE, fnm), oenv, "Dihedral Autocorrelation Function", nframes, nangles, dih, dt, mode, bAverCorr); } } /* Determine the non-zero part of the distribution */ for (first = 0; (first < maxangstat-1) && (angstat[first+1] == 0); first++) { ; } for (last = maxangstat-1; (last > 0) && (angstat[last-1] == 0); last--) { ; } aver = aver2 = 0; for (i = 0; (i < nframes); i++) { aver += RAD2DEG*aver_angle[i]; aver2 += sqr(RAD2DEG*aver_angle[i]); } aver /= (real) nframes; aver2 /= (real) nframes; aversig = sqrt(aver2-sqr(aver)); printf("Found points in the range from %d to %d (max %d)\n", first, last, maxangstat); printf(" < angle > = %g\n", aver); printf("< angle^2 > = %g\n", aver2); printf("Std. Dev. = %g\n", aversig); if (mult == 3) { sprintf(title, "Angle Distribution: %s", grpname); } else { sprintf(title, "Dihedral Distribution: %s", grpname); calc_distribution_props(maxangstat, angstat, -180.0, 0, NULL, &S2); fprintf(stderr, "Order parameter S^2 = %g\n", S2); } bPeriodic = (mult == 4) && (first == 0) && (last == maxangstat-1); out = xvgropen(opt2fn("-od", NFILE, fnm), title, "Degrees", "", oenv); if (output_env_get_print_xvgr_codes(oenv)) { fprintf(out, "@ subtitle \"average angle: %g\\So\\N\"\n", aver); } norm_fac = 1.0/(nangles*nframes*binwidth); if (bPeriodic) { maxstat = 0; for (i = first; (i <= last); i++) { maxstat = max(maxstat, angstat[i]*norm_fac); } fprintf(out, "@with g0\n"); fprintf(out, "@ world xmin -180\n"); fprintf(out, "@ world xmax 180\n"); fprintf(out, "@ world ymin 0\n"); fprintf(out, "@ world ymax %g\n", maxstat*1.05); fprintf(out, "@ xaxis tick major 60\n"); fprintf(out, "@ xaxis tick minor 30\n"); fprintf(out, "@ yaxis tick major 0.005\n"); fprintf(out, "@ yaxis tick minor 0.0025\n"); } for (i = first; (i <= last); i++) { fprintf(out, "%10g %10f\n", i*binwidth+180.0-maxang, angstat[i]*norm_fac); } if (bPeriodic) { /* print first bin again as last one */ fprintf(out, "%10g %10f\n", 180.0, angstat[0]*norm_fac); } ffclose(out); do_view(oenv, opt2fn("-od", NFILE, fnm), "-nxy"); if (bAver) { do_view(oenv, opt2fn("-ov", NFILE, fnm), "-nxy"); } thanx(stderr); return 0; }
static void histogramming(FILE *log,int nbin,gmx_residuetype_t rt, int nf,int maxchi,real **dih, int nlist,t_dlist dlist[], atom_id index[], gmx_bool bPhi,gmx_bool bPsi,gmx_bool bOmega,gmx_bool bChi, gmx_bool bNormalize,gmx_bool bSSHisto,const char *ssdump, real bfac_max,t_atoms *atoms, gmx_bool bDo_jc, const char *fn, const output_env_t oenv) { /* also gets 3J couplings and order parameters S2 */ t_karplus kkkphi[] = { { "J_NHa1", 6.51, -1.76, 1.6, -M_PI/3, 0.0, 0.0 }, { "J_NHa2", 6.51, -1.76, 1.6, M_PI/3, 0.0, 0.0 }, { "J_HaC'", 4.0, 1.1, 0.1, 0.0, 0.0, 0.0 }, { "J_NHCb", 4.7, -1.5, -0.2, M_PI/3, 0.0, 0.0 }, { "J_Ci-1Hai", 4.5, -1.3, -1.2, 2*M_PI/3, 0.0, 0.0 } }; t_karplus kkkpsi[] = { { "J_HaN", -0.88, -0.61,-0.27,M_PI/3, 0.0, 0.0 } }; t_karplus kkkchi1[] = { { "JHaHb2", 9.5, -1.6, 1.8, -M_PI/3, 0, 0.0 }, { "JHaHb3", 9.5, -1.6, 1.8, 0, 0, 0.0 } }; #define NKKKPHI asize(kkkphi) #define NKKKPSI asize(kkkpsi) #define NKKKCHI asize(kkkchi1) #define NJC (NKKKPHI+NKKKPSI+NKKKCHI) FILE *fp,*ssfp[3]={NULL,NULL,NULL}; const char *sss[3] = { "sheet", "helix", "coil" }; real S2; real *normhisto; real **Jc,**Jcsig; int ****his_aa_ss=NULL; int ***his_aa,**his_aa1,*histmp; int i,j,k,m,n,nn,Dih,nres,hindex,angle; gmx_bool bBfac,bOccup; char hisfile[256],hhisfile[256],sshisfile[256],title[256],*ss_str=NULL; char **leg; const char *residue_name; int rt_size; rt_size = gmx_residuetype_get_size(rt); if (bSSHisto) { fp = ffopen(ssdump,"r"); if(1 != fscanf(fp,"%d",&nres)) { gmx_fatal(FARGS,"Error reading from file %s",ssdump); } snew(ss_str,nres+1); if( 1 != fscanf(fp,"%s",ss_str)) { gmx_fatal(FARGS,"Error reading from file %s",ssdump); } ffclose(fp); /* Four dimensional array... Very cool */ snew(his_aa_ss,3); for(i=0; (i<3); i++) { snew(his_aa_ss[i],rt_size+1); for(j=0; (j<=rt_size); j++) { snew(his_aa_ss[i][j],edMax); for(Dih=0; (Dih<edMax); Dih++) snew(his_aa_ss[i][j][Dih],nbin+1); } } } snew(his_aa,edMax); for(Dih=0; (Dih<edMax); Dih++) { snew(his_aa[Dih],rt_size+1); for(i=0; (i<=rt_size); i++) { snew(his_aa[Dih][i],nbin+1); } } snew(histmp,nbin); snew(Jc,nlist); snew(Jcsig,nlist); for(i=0; (i<nlist); i++) { snew(Jc[i],NJC); snew(Jcsig[i],NJC); } j=0; n=0; for (Dih=0; (Dih<NONCHI+maxchi); Dih++) { for(i=0; (i<nlist); i++) { if (((Dih < edOmega) ) || ((Dih == edOmega) && (has_dihedral(edOmega,&(dlist[i])))) || ((Dih > edOmega) && (dlist[i].atm.Cn[Dih-NONCHI+3] != -1))) { make_histo(log,nf,dih[j],nbin,histmp,-M_PI,M_PI); if (bSSHisto) { /* Assume there is only one structure, the first. * Compute index in histogram. */ /* Check the atoms to see whether their B-factors are low enough * Check atoms to see their occupancy is 1. */ bBfac = bOccup = TRUE; for(nn=0; (nn<4); nn++,n++) { bBfac = bBfac && (atoms->pdbinfo[index[n]].bfac <= bfac_max); bOccup = bOccup && (atoms->pdbinfo[index[n]].occup == 1); } if (bOccup && ((bfac_max <= 0) || ((bfac_max > 0) && bBfac))) { hindex = ((dih[j][0]+M_PI)*nbin)/(2*M_PI); range_check(hindex,0,nbin); /* Assign dihedral to either of the structure determined * histograms */ switch(ss_str[dlist[i].resnr]) { case 'E': his_aa_ss[0][dlist[i].index][Dih][hindex]++; break; case 'H': his_aa_ss[1][dlist[i].index][Dih][hindex]++; break; default: his_aa_ss[2][dlist[i].index][Dih][hindex]++; break; } } else if (debug) fprintf(debug,"Res. %d has imcomplete occupancy or bfacs > %g\n", dlist[i].resnr,bfac_max); } else n += 4; switch (Dih) { case edPhi: calc_distribution_props(nbin,histmp,-M_PI,NKKKPHI,kkkphi,&S2); for(m=0; (m<NKKKPHI); m++) { Jc[i][m] = kkkphi[m].Jc; Jcsig[i][m] = kkkphi[m].Jcsig; } break; case edPsi: calc_distribution_props(nbin,histmp,-M_PI,NKKKPSI,kkkpsi,&S2); for(m=0; (m<NKKKPSI); m++) { Jc[i][NKKKPHI+m] = kkkpsi[m].Jc; Jcsig[i][NKKKPHI+m] = kkkpsi[m].Jcsig; } break; case edChi1: calc_distribution_props(nbin,histmp,-M_PI,NKKKCHI,kkkchi1,&S2); for(m=0; (m<NKKKCHI); m++) { Jc[i][NKKKPHI+NKKKPSI+m] = kkkchi1[m].Jc; Jcsig[i][NKKKPHI+NKKKPSI+m] = kkkchi1[m].Jcsig; } break; default: /* covers edOmega and higher Chis than Chi1 */ calc_distribution_props(nbin,histmp,-M_PI,0,NULL,&S2); break; } dlist[i].S2[Dih] = S2; /* Sum distribution per amino acid type as well */ for(k=0; (k<nbin); k++) { his_aa[Dih][dlist[i].index][k] += histmp[k]; histmp[k] = 0; } j++; } else { /* dihed not defined */ dlist[i].S2[Dih] = 0.0 ; } } } sfree(histmp); /* Print out Jcouplings */ fprintf(log,"\n *** J-Couplings from simulation (plus std. dev.) ***\n\n"); fprintf(log,"Residue "); for(i=0; (i<NKKKPHI); i++) fprintf(log,"%7s SD",kkkphi[i].name); for(i=0; (i<NKKKPSI); i++) fprintf(log,"%7s SD",kkkpsi[i].name); for(i=0; (i<NKKKCHI); i++) fprintf(log,"%7s SD",kkkchi1[i].name); fprintf(log,"\n"); for(i=0; (i<NJC+1); i++) fprintf(log,"------------"); fprintf(log,"\n"); for(i=0; (i<nlist); i++) { fprintf(log,"%-10s",dlist[i].name); for(j=0; (j<NJC); j++) fprintf(log," %5.2f %4.2f",Jc[i][j],Jcsig[i][j]); fprintf(log,"\n"); } fprintf(log,"\n"); /* and to -jc file... */ if (bDo_jc) { fp=xvgropen(fn,"\\S3\\NJ-Couplings from Karplus Equation","Residue", "Coupling",oenv); snew(leg,NJC); for(i=0; (i<NKKKPHI); i++){ leg[i] = strdup(kkkphi[i].name); } for(i=0; (i<NKKKPSI); i++){ leg[i+NKKKPHI]=strdup(kkkpsi[i].name); } for(i=0; (i<NKKKCHI); i++){ leg[i+NKKKPHI+NKKKPSI]=strdup(kkkchi1[i].name); } xvgr_legend(fp,NJC,(const char**)leg,oenv); fprintf(fp,"%5s ","#Res."); for(i=0; (i<NJC); i++) fprintf(fp,"%10s ",leg[i]); fprintf(fp,"\n"); for(i=0; (i<nlist); i++) { fprintf(fp,"%5d ",dlist[i].resnr); for(j=0; (j<NJC); j++) fprintf(fp," %8.3f",Jc[i][j]); fprintf(fp,"\n"); } ffclose(fp); for(i=0; (i<NJC); i++) sfree(leg[i]); } /* finished -jc stuff */ snew(normhisto,nbin); for(i=0; (i<rt_size); i++) { for(Dih=0; (Dih<edMax); Dih++){ /* First check whether something is in there */ for(j=0; (j<nbin); j++) if (his_aa[Dih][i][j] != 0) break; if ((j < nbin) && ((bPhi && (Dih==edPhi)) || (bPsi && (Dih==edPsi)) || (bOmega &&(Dih==edOmega)) || (bChi && (Dih>=edChi1)))) { if (bNormalize) normalize_histo(nbin,his_aa[Dih][i],(360.0/nbin),normhisto); residue_name = gmx_residuetype_get_name(rt,i); switch (Dih) { case edPhi: sprintf(hisfile,"histo-phi%s",residue_name); sprintf(title,"\\xf\\f{} Distribution for %s",residue_name); break; case edPsi: sprintf(hisfile,"histo-psi%s",residue_name); sprintf(title,"\\xy\\f{} Distribution for %s",residue_name); break; case edOmega: sprintf(hisfile,"histo-omega%s",residue_name); sprintf(title,"\\xw\\f{} Distribution for %s",residue_name); break; default: sprintf(hisfile,"histo-chi%d%s",Dih-NONCHI+1,residue_name); sprintf(title,"\\xc\\f{}\\s%d\\N Distribution for %s", Dih-NONCHI+1,residue_name); } strcpy(hhisfile,hisfile); strcat(hhisfile,".xvg"); fp=xvgropen(hhisfile,title,"Degrees","",oenv); fprintf(fp,"@ with g0\n"); xvgr_world(fp,-180,0,180,0.1,oenv); fprintf(fp,"# this effort to set graph size fails unless you run with -autoscale none or -autoscale y flags\n"); fprintf(fp,"@ xaxis tick on\n"); fprintf(fp,"@ xaxis tick major 90\n"); fprintf(fp,"@ xaxis tick minor 30\n"); fprintf(fp,"@ xaxis ticklabel prec 0\n"); fprintf(fp,"@ yaxis tick off\n"); fprintf(fp,"@ yaxis ticklabel off\n"); fprintf(fp,"@ type xy\n"); if (bSSHisto) { for(k=0; (k<3); k++) { sprintf(sshisfile,"%s-%s.xvg",hisfile,sss[k]); ssfp[k] = ffopen(sshisfile,"w"); } } for(j=0; (j<nbin); j++) { angle = -180 + (360/nbin)*j ; if (bNormalize) fprintf(fp,"%5d %10g\n",angle,normhisto[j]); else fprintf(fp,"%5d %10d\n",angle,his_aa[Dih][i][j]); if (bSSHisto) for(k=0; (k<3); k++) fprintf(ssfp[k],"%5d %10d\n",angle, his_aa_ss[k][i][Dih][j]); } fprintf(fp,"&\n"); ffclose(fp); if (bSSHisto) { for(k=0; (k<3); k++) { fprintf(ssfp[k],"&\n"); ffclose(ssfp[k]); } } } } } sfree(normhisto); if (bSSHisto) { /* Four dimensional array... Very cool */ for(i=0; (i<3); i++) { for(j=0; (j<=rt_size); j++) { for(Dih=0; (Dih<edMax); Dih++) sfree(his_aa_ss[i][j][Dih]); sfree(his_aa_ss[i][j]); } sfree(his_aa_ss[i]); } sfree(his_aa_ss); sfree(ss_str); } }