521ed5b782a57d410a704794e458319b53409a35
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_helix.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  * Copyright (c) 2013, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #ifdef HAVE_CONFIG_H
38 #include <config.h>
39 #endif
40
41 #include <math.h>
42
43 #include "gromacs/fileio/confio.h"
44 #include "gmx_fatal.h"
45 #include "fitahx.h"
46 #include "gromacs/fileio/futil.h"
47 #include "gstat.h"
48 #include "gromacs/fileio/g87io.h"
49 #include "hxprops.h"
50 #include "macros.h"
51 #include "maths.h"
52 #include "pbc.h"
53 #include "gromacs/fileio/tpxio.h"
54 #include "gromacs/fileio/trxio.h"
55 #include "physics.h"
56 #include "index.h"
57 #include "smalloc.h"
58 #include "gromacs/commandline/pargs.h"
59 #include <string.h>
60 #include "sysstuff.h"
61 #include "txtdump.h"
62 #include "typedefs.h"
63 #include "vec.h"
64 #include "xvgr.h"
65 #include "gmx_ana.h"
66
67 void dump_otrj(FILE *otrj, int natoms, atom_id all_index[], rvec x[],
68                real fac, rvec xav[])
69 {
70     static FILE *fp   = NULL;
71     static real  fac0 = 1.0;
72
73     int          i, ai, m;
74     real         xa, xm;
75
76     if (fp == NULL)
77     {
78         fp   = ffopen("WEDGAMMA10.DAT", "w");
79         fac0 = fac;
80     }
81     fac /= fac0;
82
83     fprintf(fp, "%10g\n", fac);
84
85     for (i = 0; (i < natoms); i++)
86     {
87         ai = all_index[i];
88         for (m = 0; (m < DIM); m++)
89         {
90             xa         = xav[ai][m];
91             xm         = x[ai][m]*fac;
92             xav[ai][m] = xa+xm;
93             x[ai][m]   = xm;
94         }
95     }
96     write_g87_ndx(otrj, natoms, all_index, x, NULL);
97 }
98
99 int gmx_helix(int argc, char *argv[])
100 {
101     const char        *desc[] = {
102         "[THISMODULE] computes all kinds of helix properties. First, the peptide",
103         "is checked to find the longest helical part, as determined by",
104         "hydrogen bonds and [GRK]phi[grk]/[GRK]psi[grk] angles.",
105         "That bit is fitted",
106         "to an ideal helix around the [IT]z[it]-axis and centered around the origin.",
107         "Then the following properties are computed:[PAR]",
108         "[BB]1.[bb] Helix radius (file [TT]radius.xvg[tt]). This is merely the",
109         "RMS deviation in two dimensions for all C[GRK]alpha[grk] atoms.",
110         "it is calculated as [SQRT]([SUM][sum][SUB]i[sub] (x^2(i)+y^2(i)))/N[sqrt] where N is the number",
111         "of backbone atoms. For an ideal helix the radius is 0.23 nm[BR]",
112         "[BB]2.[bb] Twist (file [TT]twist.xvg[tt]). The average helical angle per",
113         "residue is calculated. For an [GRK]alpha[grk]-helix it is 100 degrees,",
114         "for 3-10 helices it will be smaller, and ",
115         "for 5-helices it will be larger.[BR]",
116         "[BB]3.[bb] Rise per residue (file [TT]rise.xvg[tt]). The helical rise per",
117         "residue is plotted as the difference in [IT]z[it]-coordinate between C[GRK]alpha[grk]",
118         "atoms. For an ideal helix, this is 0.15 nm[BR]",
119         "[BB]4.[bb] Total helix length (file [TT]len-ahx.xvg[tt]). The total length",
120         "of the",
121         "helix in nm. This is simply the average rise (see above) times the",
122         "number of helical residues (see below).[BR]",
123         "[BB]5.[bb] Helix dipole, backbone only (file [TT]dip-ahx.xvg[tt]).[BR]",
124         "[BB]6.[bb] RMS deviation from ideal helix, calculated for the C[GRK]alpha[grk]",
125         "atoms only (file [TT]rms-ahx.xvg[tt]).[BR]",
126         "[BB]7.[bb] Average C[GRK]alpha[grk] - C[GRK]alpha[grk] dihedral angle (file [TT]phi-ahx.xvg[tt]).[BR]",
127         "[BB]8.[bb] Average [GRK]phi[grk] and [GRK]psi[grk] angles (file [TT]phipsi.xvg[tt]).[BR]",
128         "[BB]9.[bb] Ellipticity at 222 nm according to Hirst and Brooks.",
129         "[PAR]"
130     };
131     static const char *ppp[efhNR+2] = {
132         NULL, "RAD", "TWIST", "RISE", "LEN", "NHX", "DIP", "RMS", "CPHI",
133         "RMSA", "PHI", "PSI", "HB3", "HB4", "HB5", "CD222", NULL
134     };
135     static gmx_bool    bCheck = FALSE, bFit = TRUE, bDBG = FALSE, bEV = FALSE;
136     static int         rStart = 0, rEnd = 0, r0 = 1;
137     t_pargs            pa []  = {
138         { "-r0", FALSE, etINT, {&r0},
139           "The first residue number in the sequence" },
140         { "-q",  FALSE, etBOOL, {&bCheck},
141           "Check at every step which part of the sequence is helical" },
142         { "-F",  FALSE, etBOOL, {&bFit},
143           "Toggle fit to a perfect helix" },
144         { "-db", FALSE, etBOOL, {&bDBG},
145           "Print debug info" },
146         { "-prop", FALSE, etENUM, {ppp},
147           "Select property to weight eigenvectors with. WARNING experimental stuff" },
148         { "-ev", FALSE, etBOOL, {&bEV},
149           "Write a new 'trajectory' file for ED" },
150         { "-ahxstart", FALSE, etINT, {&rStart},
151           "First residue in helix" },
152         { "-ahxend", FALSE, etINT, {&rEnd},
153           "Last residue in helix" }
154     };
155
156     typedef struct {
157         FILE       *fp, *fp2;
158         gmx_bool    bfp2;
159         const char *filenm;
160         const char *title;
161         const char *xaxis;
162         const char *yaxis;
163         real        val;
164     } t_xvgrfile;
165
166     t_xvgrfile     xf[efhNR] = {
167         { NULL, NULL, TRUE,  "radius",  "Helix radius",               NULL, "r (nm)", 0.0 },
168         { NULL, NULL, TRUE,  "twist",   "Twist per residue",          NULL, "Angle (deg)", 0.0 },
169         { NULL, NULL, TRUE,  "rise",    "Rise per residue",           NULL, "Rise (nm)", 0.0 },
170         { NULL, NULL, FALSE, "len-ahx", "Length of the Helix",        NULL, "Length (nm)", 0.0 },
171         { NULL, NULL, FALSE, "dip-ahx", "Helix Backbone Dipole",      NULL, "rq (nm e)", 0.0 },
172         { NULL, NULL, TRUE,  "rms-ahx", "RMS Deviation from Ideal Helix", NULL, "RMS (nm)", 0.0 },
173         { NULL, NULL, FALSE, "rmsa-ahx", "Average RMSD per Residue",   "Residue", "RMS (nm)", 0.0 },
174         { NULL, NULL, FALSE,  "cd222",   "Ellipticity at 222 nm", NULL, "nm", 0.0 },
175         { NULL, NULL, TRUE,  "pprms",   "RMS Distance from \\8a\\4-helix", NULL, "deg", 0.0 },
176         { NULL, NULL, TRUE,  "caphi",   "Average Ca-Ca Dihedral",     NULL, "\\8F\\4(deg)", 0.0 },
177         { NULL, NULL, TRUE,  "phi",     "Average \\8F\\4 angles", NULL, "deg", 0.0 },
178         { NULL, NULL, TRUE,  "psi",     "Average \\8Y\\4 angles", NULL, "deg", 0.0 },
179         { NULL, NULL, TRUE,  "hb3",     "Average n-n+3 hbond length", NULL, "nm", 0.0 },
180         { NULL, NULL, TRUE,  "hb4",     "Average n-n+4 hbond length", NULL, "nm", 0.0 },
181         { NULL, NULL, TRUE,  "hb5",     "Average n-n+5 hbond length", NULL, "nm", 0.0 },
182         { NULL, NULL, FALSE,  "JCaHa",   "J-Coupling Values",        "Residue", "Hz", 0.0 },
183         { NULL, NULL, FALSE,  "helicity", "Helicity per Residue",     "Residue", "% of time", 0.0 }
184     };
185
186     output_env_t   oenv;
187     FILE          *otrj;
188     char           buf[54], prop[256];
189     t_trxstatus   *status;
190     int            natoms, nre, nres;
191     t_bb          *bb;
192     int            i, j, ai, m, nall, nbb, nca, teller, nSel = 0;
193     atom_id       *bbindex, *caindex, *allindex;
194     t_topology    *top;
195     int            ePBC;
196     rvec          *x, *xref, *xav;
197     real           t;
198     real           rms, fac;
199     matrix         box;
200     gmx_rmpbc_t    gpbc = NULL;
201     gmx_bool       bRange;
202     t_filenm       fnm[] = {
203         { efTPX, NULL,  NULL,   ffREAD  },
204         { efNDX, NULL,  NULL,   ffREAD  },
205         { efTRX, "-f",  NULL,   ffREAD  },
206         { efG87, "-to", NULL,   ffOPTWR },
207         { efSTO, "-cz", "zconf", ffWRITE },
208         { efSTO, "-co", "waver", ffWRITE }
209     };
210 #define NFILE asize(fnm)
211
212     if (!parse_common_args(&argc, argv, PCA_CAN_VIEW | PCA_CAN_TIME | PCA_BE_NICE,
213                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
214     {
215         return 0;
216     }
217
218     bRange = (opt2parg_bSet("-ahxstart", asize(pa), pa) &&
219               opt2parg_bSet("-ahxend", asize(pa), pa));
220
221     top = read_top(ftp2fn(efTPX, NFILE, fnm), &ePBC);
222
223     natoms = read_first_x(oenv, &status, opt2fn("-f", NFILE, fnm), &t, &x, box);
224
225     if (opt2bSet("-to", NFILE, fnm))
226     {
227         otrj = opt2FILE("-to", NFILE, fnm, "w");
228         strcpy(prop, ppp[0]);
229         fprintf(otrj, "%s Weighted Trajectory: %d atoms, NO box\n", prop, natoms);
230     }
231     else
232     {
233         otrj = NULL;
234     }
235
236     if (natoms != top->atoms.nr)
237     {
238         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)",
239                   top->atoms.nr, natoms);
240     }
241
242     bb = mkbbind(ftp2fn(efNDX, NFILE, fnm), &nres, &nbb, r0, &nall, &allindex,
243                  top->atoms.atomname, top->atoms.atom, top->atoms.resinfo);
244     snew(bbindex, natoms);
245     snew(caindex, nres);
246
247     fprintf(stderr, "nall=%d\n", nall);
248
249     /* Open output files, default x-axis is time */
250     for (i = 0; (i < efhNR); i++)
251     {
252         sprintf(buf, "%s.xvg", xf[i].filenm);
253         remove(buf);
254         xf[i].fp = xvgropen(buf, xf[i].title,
255                             xf[i].xaxis ? xf[i].xaxis : "Time (ps)",
256                             xf[i].yaxis, oenv);
257         if (xf[i].bfp2)
258         {
259             sprintf(buf, "%s.out", xf[i].filenm);
260             remove(buf);
261             xf[i].fp2 = ffopen(buf, "w");
262         }
263     }
264
265     /* Read reference frame from tpx file to compute helix length */
266     snew(xref, top->atoms.nr);
267     read_tpx(ftp2fn(efTPX, NFILE, fnm),
268              NULL, NULL, &natoms, xref, NULL, NULL, NULL);
269     calc_hxprops(nres, bb, xref);
270     do_start_end(nres, bb, &nbb, bbindex, &nca, caindex, bRange, rStart, rEnd);
271     sfree(xref);
272     if (bDBG)
273     {
274         fprintf(stderr, "nca=%d, nbb=%d\n", nca, nbb);
275         pr_bb(stdout, nres, bb);
276     }
277
278     gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms);
279
280     snew(xav, natoms);
281     teller = 0;
282     do
283     {
284         if ((teller++ % 10) == 0)
285         {
286             fprintf(stderr, "\rt=%.2f", t);
287         }
288         gmx_rmpbc(gpbc, natoms, box, x);
289
290
291         calc_hxprops(nres, bb, x);
292         if (bCheck)
293         {
294             do_start_end(nres, bb, &nbb, bbindex, &nca, caindex, FALSE, 0, 0);
295         }
296
297         if (nca >= 5)
298         {
299             rms = fit_ahx(nres, bb, natoms, nall, allindex, x, nca, caindex, bFit);
300
301             if (teller == 1)
302             {
303                 write_sto_conf(opt2fn("-cz", NFILE, fnm), "Helix fitted to Z-Axis",
304                                &(top->atoms), x, NULL, ePBC, box);
305             }
306
307             xf[efhRAD].val   = radius(xf[efhRAD].fp2, nca, caindex, x);
308             xf[efhTWIST].val = twist(nca, caindex, x);
309             xf[efhRISE].val  = rise(nca, caindex, x);
310             xf[efhLEN].val   = ahx_len(nca, caindex, x);
311             xf[efhCD222].val = ellipticity(nres, bb);
312             xf[efhDIP].val   = dip(nbb, bbindex, x, top->atoms.atom);
313             xf[efhRMS].val   = rms;
314             xf[efhCPHI].val  = ca_phi(nca, caindex, x);
315             xf[efhPPRMS].val = pprms(xf[efhPPRMS].fp2, nres, bb);
316
317             for (j = 0; (j <= efhCPHI); j++)
318             {
319                 fprintf(xf[j].fp,   "%10g  %10g\n", t, xf[j].val);
320             }
321
322             av_phipsi(xf[efhPHI].fp, xf[efhPSI].fp, xf[efhPHI].fp2, xf[efhPSI].fp2,
323                       t, nres, bb);
324             av_hblen(xf[efhHB3].fp, xf[efhHB3].fp2,
325                      xf[efhHB4].fp, xf[efhHB4].fp2,
326                      xf[efhHB5].fp, xf[efhHB5].fp2,
327                      t, nres, bb);
328
329             if (otrj)
330             {
331                 dump_otrj(otrj, nall, allindex, x, xf[nSel].val, xav);
332             }
333         }
334     }
335     while (read_next_x(oenv, status, &t, x, box));
336     fprintf(stderr, "\n");
337
338     gmx_rmpbc_done(gpbc);
339
340     close_trj(status);
341
342     if (otrj)
343     {
344         ffclose(otrj);
345         fac = 1.0/teller;
346         for (i = 0; (i < nall); i++)
347         {
348             ai = allindex[i];
349             for (m = 0; (m < DIM); m++)
350             {
351                 xav[ai][m] *= fac;
352             }
353         }
354         write_sto_conf_indexed(opt2fn("-co", NFILE, fnm),
355                                "Weighted and Averaged conformation",
356                                &(top->atoms), xav, NULL, ePBC, box, nall, allindex);
357     }
358
359     for (i = 0; (i < nres); i++)
360     {
361         if (bb[i].nrms > 0)
362         {
363             fprintf(xf[efhRMSA].fp, "%10d  %10g\n", r0+i, bb[i].rmsa/bb[i].nrms);
364         }
365         fprintf(xf[efhAHX].fp, "%10d  %10g\n", r0+i, (bb[i].nhx*100.0)/(real )teller);
366         fprintf(xf[efhJCA].fp, "%10d  %10g\n",
367                 r0+i, 140.3+(bb[i].jcaha/(double)teller));
368     }
369
370     for (i = 0; (i < efhNR); i++)
371     {
372         ffclose(xf[i].fp);
373         if (xf[i].bfp2)
374         {
375             ffclose(xf[i].fp2);
376         }
377         do_view(oenv, xf[i].filenm, "-nxy");
378     }
379
380     return 0;
381 }