905fa4e1b45320b647fb1e969e4a5ce83eda3cd2
[alexxy/gromacs.git] / src / tools / gmx_helix.c
1 /*
2  * 
3  *                This source code is part of
4  * 
5  *                 G   R   O   M   A   C   S
6  * 
7  *          GROningen MAchine for Chemical Simulations
8  * 
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  * 
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  * 
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  * 
30  * For more info, check our website at http://www.gromacs.org
31  * 
32  * And Hey:
33  * Green Red Orange Magenta Azure Cyan Skyblue
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <math.h>
40
41 #include "confio.h"
42 #include "copyrite.h"
43 #include "gmx_fatal.h"
44 #include "fitahx.h"
45 #include "futil.h"
46 #include "gstat.h"
47 #include "wgms.h"
48 #include "hxprops.h"
49 #include "macros.h"
50 #include "maths.h"
51 #include "pbc.h"
52 #include "tpxio.h"
53 #include "physics.h"
54 #include "index.h"
55 #include "smalloc.h"
56 #include "statutil.h"
57 #include <string.h>
58 #include "sysstuff.h"
59 #include "txtdump.h"
60 #include "typedefs.h"
61 #include "vec.h"
62 #include "xvgr.h"
63 #include "gmx_ana.h"
64
65 void dump_otrj(FILE *otrj,int natoms,atom_id all_index[],rvec x[],
66                real fac,rvec xav[])
67 {
68   static FILE *fp=NULL;
69   static real fac0=1.0;
70   
71   int  i,ai,m;
72   real xa,xm;
73   
74   if (fp==NULL) {
75     fp=ffopen("WEDGAMMA10.DAT","w");
76     fac0=fac;
77   }
78   fac/=fac0;
79   
80   fprintf(fp,"%10g\n",fac);
81   
82   for(i=0; (i<natoms); i++) {
83     ai=all_index[i];
84     for(m=0; (m<DIM); m++) {
85       xa = xav[ai][m];
86       xm = x[ai][m]*fac;
87       xav[ai][m] = xa+xm;
88       x[ai][m]   = xm;
89     }
90   }
91   write_gms_ndx(otrj,natoms,all_index,x,NULL);
92 }
93
94 int gmx_helix(int argc,char *argv[])
95 {
96   const char *desc[] = {
97     "[TT]g_helix[tt] computes all kinds of helix properties. First, the peptide",
98     "is checked to find the longest helical part, as determined by",
99     "hydrogen bonds and [GRK]phi[grk]/[GRK]psi[grk] angles.",
100     "That bit is fitted",
101     "to an ideal helix around the [IT]z[it]-axis and centered around the origin.",
102     "Then the following properties are computed:[PAR]",
103     "[BB]1.[bb] Helix radius (file [TT]radius.xvg[tt]). This is merely the",
104     "RMS deviation in two dimensions for all C[GRK]alpha[grk] atoms.",
105     "it is calculated as [SQRT]([SUM][sum][SUB]i[sub] (x^2(i)+y^2(i)))/N[sqrt] where N is the number",
106     "of backbone atoms. For an ideal helix the radius is 0.23 nm[BR]",
107     "[BB]2.[bb] Twist (file [TT]twist.xvg[tt]). The average helical angle per",
108     "residue is calculated. For an [GRK]alpha[grk]-helix it is 100 degrees,",
109     "for 3-10 helices it will be smaller, and ", 
110     "for 5-helices it will be larger.[BR]",
111     "[BB]3.[bb] Rise per residue (file [TT]rise.xvg[tt]). The helical rise per", 
112     "residue is plotted as the difference in [IT]z[it]-coordinate between C[GRK]alpha[grk]", 
113     "atoms. For an ideal helix, this is 0.15 nm[BR]",
114     "[BB]4.[bb] Total helix length (file [TT]len-ahx.xvg[tt]). The total length", 
115     "of the", 
116     "helix in nm. This is simply the average rise (see above) times the",  
117     "number of helical residues (see below).[BR]",
118     "[BB]5.[bb] Number of helical residues (file [TT]n-ahx.xvg[tt]). The title says",
119     "it all.[BR]",
120     "[BB]6.[bb] Helix dipole, backbone only (file [TT]dip-ahx.xvg[tt]).[BR]",
121     "[BB]7.[bb] RMS deviation from ideal helix, calculated for the C[GRK]alpha[grk]",
122     "atoms only (file [TT]rms-ahx.xvg[tt]).[BR]",
123     "[BB]8.[bb] Average C[GRK]alpha[grk] - C[GRK]alpha[grk] dihedral angle (file [TT]phi-ahx.xvg[tt]).[BR]",
124     "[BB]9.[bb] Average [GRK]phi[grk] and [GRK]psi[grk] angles (file [TT]phipsi.xvg[tt]).[BR]",
125     "[BB]10.[bb] Ellipticity at 222 nm according to Hirst and Brooks.",
126     "[PAR]"
127   };
128   static const char *ppp[efhNR+2] = { 
129     NULL, "RAD", "TWIST", "RISE", "LEN", "NHX", "DIP", "RMS", "CPHI", 
130     "RMSA", "PHI", "PSI", "HB3", "HB4", "HB5", "CD222", NULL
131   };
132   static gmx_bool bCheck=FALSE,bFit=TRUE,bDBG=FALSE,bEV=FALSE;
133   static int  rStart=0,rEnd=0,r0=1;
134   t_pargs pa [] = {
135     { "-r0", FALSE, etINT, {&r0},
136       "The first residue number in the sequence" },
137     { "-q",  FALSE, etBOOL,{&bCheck},
138       "Check at every step which part of the sequence is helical" },
139     { "-F",  FALSE, etBOOL,{&bFit},
140       "Toggle fit to a perfect helix" },
141     { "-db", FALSE, etBOOL,{&bDBG},
142       "Print debug info" },
143     { "-prop", FALSE, etENUM, {ppp},
144       "Select property to weight eigenvectors with. WARNING experimental stuff" },
145     { "-ev", FALSE, etBOOL,{&bEV},
146       "Write a new 'trajectory' file for ED" },
147     { "-ahxstart", FALSE, etINT, {&rStart},
148       "First residue in helix" },
149     { "-ahxend", FALSE, etINT, {&rEnd},
150       "Last residue in helix" }
151   };
152
153   typedef struct {
154     FILE *fp,*fp2;
155     gmx_bool bfp2;
156     const char *filenm;
157     const char *title;
158     const char *xaxis;
159     const char *yaxis;
160     real val;
161   } t_xvgrfile;
162   
163   t_xvgrfile xf[efhNR] = {
164     { NULL, NULL, TRUE,  "radius",  "Helix radius",               NULL, "r (nm)" , 0.0 },
165     { NULL, NULL, TRUE,  "twist",   "Twist per residue",          NULL, "Angle (deg)", 0.0 },
166     { NULL, NULL, TRUE,  "rise",    "Rise per residue",           NULL, "Rise (nm)", 0.0 },
167     { NULL, NULL, FALSE, "len-ahx", "Length of the Helix",        NULL, "Length (nm)", 0.0 },
168     { NULL, NULL, FALSE, "dip-ahx", "Helix Backbone Dipole",      NULL, "rq (nm e)", 0.0 },
169     { NULL, NULL, TRUE,  "rms-ahx", "RMS Deviation from Ideal Helix", NULL, "RMS (nm)", 0.0 },
170     { NULL, NULL, FALSE, "rmsa-ahx","Average RMSD per Residue",   "Residue", "RMS (nm)", 0.0 },
171     { NULL, NULL,FALSE,  "cd222",   "Ellipticity at 222 nm", NULL, "nm", 0.0 },
172     { NULL, NULL, TRUE,  "pprms",   "RMS Distance from \\8a\\4-helix", NULL, "deg" , 0.0 },
173     { NULL, NULL, TRUE,  "caphi",   "Average Ca-Ca Dihedral",     NULL, "\\8F\\4(deg)", 0.0 },
174     { NULL, NULL, TRUE,  "phi",     "Average \\8F\\4 angles", NULL, "deg" , 0.0 },
175     { NULL, NULL, TRUE,  "psi",     "Average \\8Y\\4 angles", NULL, "deg" , 0.0 },
176     { NULL, NULL, TRUE,  "hb3",     "Average n-n+3 hbond length", NULL, "nm" , 0.0 },
177     { NULL, NULL, TRUE,  "hb4",     "Average n-n+4 hbond length", NULL, "nm" , 0.0 },
178     { NULL, NULL, TRUE,  "hb5",     "Average n-n+5 hbond length", NULL, "nm" , 0.0 },
179     { NULL, NULL,FALSE,  "JCaHa",   "J-Coupling Values",        "Residue", "Hz" , 0.0 },
180     { NULL, NULL,FALSE,  "helicity","Helicity per Residue",     "Residue", "% of time" , 0.0 }
181   };
182  
183   output_env_t oenv;
184   FILE       *otrj;
185   char       buf[54],prop[256];
186   t_trxstatus *status;
187   int        natoms,nre,nres;
188   t_bb       *bb;
189   int        i,j,ai,m,nall,nbb,nca,teller,nSel=0;
190   atom_id    *bbindex,*caindex,*allindex;
191   t_topology *top;
192   int        ePBC;
193   rvec       *x,*xref,*xav;
194   real       t;
195   real       rms,fac;
196   matrix     box;
197   gmx_rmpbc_t  gpbc=NULL;
198   gmx_bool       bRange;
199   t_filenm  fnm[] = {
200     { efTPX, NULL,  NULL,   ffREAD  },
201     { efNDX, NULL,  NULL,   ffREAD  },
202     { efTRX, "-f",  NULL,   ffREAD  },
203     { efG87, "-to", NULL,   ffOPTWR },
204     { efSTO, "-cz", "zconf",ffWRITE },
205     { efSTO, "-co", "waver",ffWRITE }
206   };
207 #define NFILE asize(fnm)
208
209   CopyRight(stderr,argv[0]);
210   parse_common_args(&argc,argv,PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
211                     NFILE,fnm,asize(pa),pa,asize(desc),desc,0,NULL,&oenv);
212   
213   bRange=(opt2parg_bSet("-ahxstart",asize(pa),pa) &&
214           opt2parg_bSet("-ahxend",asize(pa),pa));
215                         
216   top=read_top(ftp2fn(efTPX,NFILE,fnm),&ePBC);
217   
218   natoms=read_first_x(oenv,&status,opt2fn("-f",NFILE,fnm),&t,&x,box);
219
220   if (opt2bSet("-to",NFILE,fnm)) {
221     otrj=opt2FILE("-to",NFILE,fnm,"w");
222     strcpy(prop,ppp[0]);
223     fprintf(otrj,"%s Weighted Trajectory: %d atoms, NO box\n",prop,natoms);
224   }
225   else
226     otrj=NULL;
227     
228   if (natoms != top->atoms.nr)
229     gmx_fatal(FARGS,"Sorry can only run when the number of atoms in the run input file (%d) is equal to the number in the trajectory (%d)",
230             top->atoms.nr,natoms);
231             
232   bb=mkbbind(ftp2fn(efNDX,NFILE,fnm),&nres,&nbb,r0,&nall,&allindex,
233              top->atoms.atomname,top->atoms.atom,top->atoms.resinfo);
234   snew(bbindex,natoms);
235   snew(caindex,nres);
236   
237   fprintf(stderr,"nall=%d\n",nall);
238     
239   /* Open output files, default x-axis is time */
240   for(i=0; (i<efhNR); i++) {
241     sprintf(buf,"%s.xvg",xf[i].filenm);
242     remove(buf);
243     xf[i].fp=xvgropen(buf,xf[i].title,
244                       xf[i].xaxis ? xf[i].xaxis : "Time (ps)",
245                       xf[i].yaxis,oenv);
246     if (xf[i].bfp2) {
247       sprintf(buf,"%s.out",xf[i].filenm);
248       remove(buf);
249       xf[i].fp2=ffopen(buf,"w");
250     }
251   }
252
253   /* Read reference frame from tpx file to compute helix length */
254   snew(xref,top->atoms.nr);
255   read_tpx(ftp2fn(efTPX,NFILE,fnm),
256            NULL,NULL,&natoms,xref,NULL,NULL,NULL);
257   calc_hxprops(nres,bb,xref,box);
258   do_start_end(nres,bb,xref,&nbb,bbindex,&nca,caindex,bRange,rStart,rEnd);
259   sfree(xref);
260   if (bDBG) {
261     fprintf(stderr,"nca=%d, nbb=%d\n",nca,nbb);
262     pr_bb(stdout,nres,bb);
263   }
264   
265   gpbc = gmx_rmpbc_init(&top->idef,ePBC,natoms,box);
266
267   snew(xav,natoms);
268   teller=0;
269   do {
270     if ((teller++ % 10) == 0)
271       fprintf(stderr,"\rt=%.2f",t);
272     gmx_rmpbc(gpbc,natoms,box,x);
273
274     
275     calc_hxprops(nres,bb,x,box);
276     if (bCheck)
277       do_start_end(nres,bb,x,&nbb,bbindex,&nca,caindex,FALSE,0,0);
278     
279     if (nca >= 5) {
280       rms=fit_ahx(nres,bb,natoms,nall,allindex,x,nca,caindex,box,bFit);
281       
282       if (teller == 1) {
283         write_sto_conf(opt2fn("-cz",NFILE,fnm),"Helix fitted to Z-Axis",
284                        &(top->atoms),x,NULL,ePBC,box);
285       }
286             
287       xf[efhRAD].val   = radius(xf[efhRAD].fp2,nca,caindex,x);
288       xf[efhTWIST].val = twist(xf[efhTWIST].fp2,nca,caindex,x);
289       xf[efhRISE].val  = rise(nca,caindex,x);
290       xf[efhLEN].val   = ahx_len(nca,caindex,x,box);
291       xf[efhCD222].val = ellipticity(nres,bb);
292       xf[efhDIP].val   = dip(nbb,bbindex,x,top->atoms.atom);
293       xf[efhRMS].val   = rms;
294       xf[efhCPHI].val  = ca_phi(nca,caindex,x,box);
295       xf[efhPPRMS].val = pprms(xf[efhPPRMS].fp2,nres,bb);
296       
297       for(j=0; (j<=efhCPHI); j++)
298         fprintf(xf[j].fp,   "%10g  %10g\n",t,xf[j].val);
299       
300       av_phipsi(xf[efhPHI].fp,xf[efhPSI].fp,xf[efhPHI].fp2,xf[efhPSI].fp2,
301                 t,nres,bb);
302       av_hblen(xf[efhHB3].fp,xf[efhHB3].fp2,
303                xf[efhHB4].fp,xf[efhHB4].fp2,
304                xf[efhHB5].fp,xf[efhHB5].fp2,
305                t,nres,bb);
306       
307       if (otrj) 
308         dump_otrj(otrj,nall,allindex,x,xf[nSel].val,xav);
309     }
310   } while (read_next_x(oenv,status,&t,natoms,x,box));
311   fprintf(stderr,"\n");
312   
313   gmx_rmpbc_done(gpbc);
314
315   close_trj(status);
316
317   if (otrj) {
318     ffclose(otrj);
319     fac=1.0/teller;
320     for(i=0; (i<nall); i++) {
321       ai=allindex[i];
322       for(m=0; (m<DIM); m++)
323         xav[ai][m]*=fac;
324     }
325     write_sto_conf_indexed(opt2fn("-co",NFILE,fnm),
326                            "Weighted and Averaged conformation",
327                            &(top->atoms),xav,NULL,ePBC,box,nall,allindex);
328   }
329   
330   for(i=0; (i<nres); i++) {
331     if (bb[i].nrms > 0) {
332       fprintf(xf[efhRMSA].fp,"%10d  %10g\n",r0+i,bb[i].rmsa/bb[i].nrms);
333     }
334     fprintf(xf[efhAHX].fp,"%10d  %10g\n",r0+i,(bb[i].nhx*100.0)/(real )teller);
335     fprintf(xf[efhJCA].fp,"%10d  %10g\n",
336             r0+i,140.3+(bb[i].jcaha/(double)teller));
337   }
338   
339   for(i=0; (i<efhNR); i++) {
340     ffclose(xf[i].fp);
341     if (xf[i].bfp2)
342       ffclose(xf[i].fp2);
343     do_view(oenv,xf[i].filenm,"-nxy");
344   }
345   
346   thanx(stderr);
347   
348   return 0;
349 }
350