9cb602d93a71dd7f5a1659b362c29dce9bc4526f
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_bundle.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,2014, 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 #include "gmxpre.h"
38
39 #include "config.h"
40 #include <math.h>
41 #include <string.h>
42
43 #include "gromacs/commandline/pargs.h"
44 #include "gromacs/legacyheaders/typedefs.h"
45 #include "gromacs/utility/smalloc.h"
46 #include "gromacs/legacyheaders/macros.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/utility/futil.h"
49 #include "gromacs/topology/index.h"
50 #include "gromacs/fileio/xvgr.h"
51 #include "gromacs/pbcutil/rmpbc.h"
52 #include "gromacs/fileio/tpxio.h"
53 #include "gromacs/fileio/trxio.h"
54 #include "gromacs/math/units.h"
55 #include "gmx_ana.h"
56
57 #include "gromacs/utility/fatalerror.h"
58
59 #define MAX_ENDS 3
60
61 typedef struct {
62     int   n;
63     int   nend;
64     rvec *end[MAX_ENDS];
65     rvec *mid;
66     rvec *dir;
67     real *len;
68 } t_bundle;
69
70 static void rotate_ends(t_bundle *bun, rvec axis, int c0, int c1)
71 {
72     int  end, i;
73     rvec ax, tmp;
74
75     unitv(axis, ax);
76     for (end = 0; end < bun->nend; end++)
77     {
78         for (i = 0; i < bun->n; i++)
79         {
80             copy_rvec(bun->end[end][i], tmp);
81             bun->end[end][i][c0] = ax[c1]*tmp[c0] - ax[c0]*tmp[c1];
82             bun->end[end][i][c1] = ax[c0]*tmp[c0] + ax[c1]*tmp[c1];
83         }
84     }
85     copy_rvec(axis, tmp);
86     axis[c0] = ax[c1]*tmp[c0] - ax[c0]*tmp[c1];
87     axis[c1] = ax[c0]*tmp[c0] + ax[c1]*tmp[c1];
88 }
89
90 static void calc_axes(rvec x[], t_atom atom[], int gnx[], atom_id *index[],
91                       gmx_bool bRot, t_bundle *bun)
92 {
93     int   end, i, div, d;
94     real *mtot, m;
95     rvec  axis[MAX_ENDS], cent;
96
97     snew(mtot, bun->n);
98
99     for (end = 0; end < bun->nend; end++)
100     {
101         for (i = 0; i < bun->n; i++)
102         {
103             clear_rvec(bun->end[end][i]);
104             mtot[i] = 0;
105         }
106         div = gnx[end]/bun->n;
107         for (i = 0; i < gnx[end]; i++)
108         {
109             m = atom[index[end][i]].m;
110             for (d = 0; d < DIM; d++)
111             {
112                 bun->end[end][i/div][d] += m*x[index[end][i]][d];
113             }
114             mtot[i/div] += m;
115         }
116         clear_rvec(axis[end]);
117         for (i = 0; i < bun->n; i++)
118         {
119             svmul(1.0/mtot[i], bun->end[end][i], bun->end[end][i]);
120             rvec_inc(axis[end], bun->end[end][i]);
121         }
122         svmul(1.0/bun->n, axis[end], axis[end]);
123     }
124     sfree(mtot);
125
126     rvec_add(axis[0], axis[1], cent);
127     svmul(0.5, cent, cent);
128     /* center the bundle on the origin */
129     for (end = 0; end < bun->nend; end++)
130     {
131         rvec_dec(axis[end], cent);
132         for (i = 0; i < bun->n; i++)
133         {
134             rvec_dec(bun->end[end][i], cent);
135         }
136     }
137     if (bRot)
138     {
139         /* rotate the axis parallel to the z-axis */
140         rotate_ends(bun, axis[0], YY, ZZ);
141         rotate_ends(bun, axis[0], XX, ZZ);
142     }
143     for (i = 0; i < bun->n; i++)
144     {
145         rvec_add(bun->end[0][i], bun->end[1][i], bun->mid[i]);
146         svmul(0.5, bun->mid[i], bun->mid[i]);
147         rvec_sub(bun->end[0][i], bun->end[1][i], bun->dir[i]);
148         bun->len[i] = norm(bun->dir[i]);
149         unitv(bun->dir[i], bun->dir[i]);
150     }
151 }
152
153 static void dump_axes(t_trxstatus *status, t_trxframe *fr, t_atoms *outat,
154                       t_bundle *bun)
155 {
156     t_trxframe   frout;
157     static rvec *xout = NULL;
158     int          i, end;
159
160     if (xout == NULL)
161     {
162         snew(xout, outat->nr);
163     }
164
165     for (i = 0; i < bun->n; i++)
166     {
167         copy_rvec(bun->end[0][i], xout[3*i]);
168         if (bun->nend >= 3)
169         {
170             copy_rvec(bun->end[2][i], xout[3*i+1]);
171         }
172         else
173         {
174             copy_rvec(bun->mid[i], xout[3*i+1]);
175         }
176         copy_rvec(bun->end[1][i], xout[3*i+2]);
177     }
178     frout        = *fr;
179     frout.bV     = FALSE;
180     frout.bF     = FALSE;
181     frout.bBox   = FALSE;
182     frout.bAtoms = TRUE;
183     frout.natoms = outat->nr;
184     frout.atoms  = outat;
185     frout.x      = xout;
186     write_trxframe(status, &frout, NULL);
187 }
188
189 int gmx_bundle(int argc, char *argv[])
190 {
191     const char     *desc[] = {
192         "[THISMODULE] analyzes bundles of axes. The axes can be for instance",
193         "helix axes. The program reads two index groups and divides both",
194         "of them in [TT]-na[tt] parts. The centers of mass of these parts",
195         "define the tops and bottoms of the axes.",
196         "Several quantities are written to file:",
197         "the axis length, the distance and the z-shift of the axis mid-points",
198         "with respect to the average center of all axes, the total tilt,",
199         "the radial tilt and the lateral tilt with respect to the average axis.",
200         "[PAR]",
201         "With options [TT]-ok[tt], [TT]-okr[tt] and [TT]-okl[tt] the total,",
202         "radial and lateral kinks of the axes are plotted. An extra index",
203         "group of kink atoms is required, which is also divided into [TT]-na[tt]",
204         "parts. The kink angle is defined as the angle between the kink-top and",
205         "the bottom-kink vectors.",
206         "[PAR]",
207         "With option [TT]-oa[tt] the top, mid (or kink when [TT]-ok[tt] is set)",
208         "and bottom points of each axis",
209         "are written to a [TT].pdb[tt] file each frame. The residue numbers correspond",
210         "to the axis numbers. When viewing this file with Rasmol, use the",
211         "command line option [TT]-nmrpdb[tt], and type [TT]set axis true[tt] to",
212         "display the reference axis."
213     };
214     static int      n    = 0;
215     static gmx_bool bZ   = FALSE;
216     t_pargs         pa[] = {
217         { "-na", FALSE, etINT, {&n},
218           "Number of axes" },
219         { "-z", FALSE, etBOOL, {&bZ},
220           "Use the [IT]z[it]-axis as reference instead of the average axis" }
221     };
222     FILE           *out, *flen, *fdist, *fz, *ftilt, *ftiltr, *ftiltl;
223     FILE           *fkink = NULL, *fkinkr = NULL, *fkinkl = NULL;
224     t_trxstatus    *status;
225     t_trxstatus    *fpdb;
226     t_topology      top;
227     int             ePBC;
228     rvec           *xtop;
229     matrix          box;
230     t_trxframe      fr;
231     t_atoms         outatoms;
232     real            t, comp;
233     int             natoms;
234     char           *grpname[MAX_ENDS], title[256];
235     /* FIXME: The constness should not be cast away */
236     char           *anm = (char *)"CA", *rnm = (char *)"GLY";
237     int             i, j, gnx[MAX_ENDS];
238     atom_id        *index[MAX_ENDS];
239     t_bundle        bun;
240     gmx_bool        bKink;
241     rvec            va, vb, vc, vr, vl;
242     output_env_t    oenv;
243     gmx_rmpbc_t     gpbc = NULL;
244
245 #define NLEG asize(leg)
246     t_filenm fnm[] = {
247         { efTRX, "-f", NULL, ffREAD },
248         { efTPS, NULL, NULL, ffREAD },
249         { efNDX, NULL, NULL, ffOPTRD },
250         { efXVG, "-ol", "bun_len", ffWRITE },
251         { efXVG, "-od", "bun_dist", ffWRITE },
252         { efXVG, "-oz", "bun_z", ffWRITE },
253         { efXVG, "-ot", "bun_tilt", ffWRITE },
254         { efXVG, "-otr", "bun_tiltr", ffWRITE },
255         { efXVG, "-otl", "bun_tiltl", ffWRITE },
256         { efXVG, "-ok", "bun_kink", ffOPTWR },
257         { efXVG, "-okr", "bun_kinkr", ffOPTWR },
258         { efXVG, "-okl", "bun_kinkl", ffOPTWR },
259         { efPDB, "-oa", "axes", ffOPTWR }
260     };
261 #define NFILE asize(fnm)
262
263     if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_TIME_UNIT,
264                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
265     {
266         return 0;
267     }
268
269     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xtop, NULL, box, TRUE);
270
271     bKink = opt2bSet("-ok", NFILE, fnm) || opt2bSet("-okr", NFILE, fnm)
272         || opt2bSet("-okl", NFILE, fnm);
273     if (bKink)
274     {
275         bun.nend = 3;
276     }
277     else
278     {
279         bun.nend = 2;
280     }
281
282     fprintf(stderr, "Select a group of top and a group of bottom ");
283     if (bKink)
284     {
285         fprintf(stderr, "and a group of kink ");
286     }
287     fprintf(stderr, "atoms\n");
288     get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), bun.nend,
289               gnx, index, grpname);
290
291     if (n <= 0 || gnx[0] % n || gnx[1] % n || (bKink && gnx[2] % n))
292     {
293         gmx_fatal(FARGS,
294                   "The size of one of your index groups is not a multiple of n");
295     }
296     bun.n = n;
297     snew(bun.end[0], n);
298     snew(bun.end[1], n);
299     if (bKink)
300     {
301         snew(bun.end[2], n);
302     }
303     snew(bun.mid, n);
304     snew(bun.dir, n);
305     snew(bun.len, n);
306
307     flen   = xvgropen(opt2fn("-ol", NFILE, fnm), "Axis lengths",
308                       output_env_get_xvgr_tlabel(oenv), "(nm)", oenv);
309     fdist  = xvgropen(opt2fn("-od", NFILE, fnm), "Distance of axis centers",
310                       output_env_get_xvgr_tlabel(oenv), "(nm)", oenv);
311     fz     = xvgropen(opt2fn("-oz", NFILE, fnm), "Z-shift of axis centers",
312                       output_env_get_xvgr_tlabel(oenv), "(nm)", oenv);
313     ftilt  = xvgropen(opt2fn("-ot", NFILE, fnm), "Axis tilts",
314                       output_env_get_xvgr_tlabel(oenv), "(degrees)", oenv);
315     ftiltr = xvgropen(opt2fn("-otr", NFILE, fnm), "Radial axis tilts",
316                       output_env_get_xvgr_tlabel(oenv), "(degrees)", oenv);
317     ftiltl = xvgropen(opt2fn("-otl", NFILE, fnm), "Lateral axis tilts",
318                       output_env_get_xvgr_tlabel(oenv), "(degrees)", oenv);
319
320     if (bKink)
321     {
322         fkink  = xvgropen(opt2fn("-ok", NFILE, fnm), "Kink angles",
323                           output_env_get_xvgr_tlabel(oenv), "(degrees)", oenv);
324         fkinkr = xvgropen(opt2fn("-okr", NFILE, fnm), "Radial kink angles",
325                           output_env_get_xvgr_tlabel(oenv), "(degrees)", oenv);
326         if (output_env_get_print_xvgr_codes(oenv))
327         {
328             fprintf(fkinkr, "@ subtitle \"+ = ) (   - = ( )\"\n");
329         }
330         fkinkl = xvgropen(opt2fn("-okl", NFILE, fnm), "Lateral kink angles",
331                           output_env_get_xvgr_tlabel(oenv), "(degrees)", oenv);
332     }
333
334     if (opt2bSet("-oa", NFILE, fnm))
335     {
336         init_t_atoms(&outatoms, 3*n, FALSE);
337         outatoms.nr = 3*n;
338         for (i = 0; i < 3*n; i++)
339         {
340             outatoms.atomname[i]       = &anm;
341             outatoms.atom[i].resind    = i/3;
342             outatoms.resinfo[i/3].name = &rnm;
343             outatoms.resinfo[i/3].nr   = i/3 + 1;
344             outatoms.resinfo[i/3].ic   = ' ';
345         }
346         fpdb = open_trx(opt2fn("-oa", NFILE, fnm), "w");
347     }
348     else
349     {
350         fpdb = NULL;
351     }
352
353     read_first_frame(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &fr, TRX_NEED_X);
354     gpbc = gmx_rmpbc_init(&top.idef, ePBC, fr.natoms);
355
356     do
357     {
358         gmx_rmpbc_trxfr(gpbc, &fr);
359         calc_axes(fr.x, top.atoms.atom, gnx, index, !bZ, &bun);
360         t = output_env_conv_time(oenv, fr.time);
361         fprintf(flen, " %10g", t);
362         fprintf(fdist, " %10g", t);
363         fprintf(fz, " %10g", t);
364         fprintf(ftilt, " %10g", t);
365         fprintf(ftiltr, " %10g", t);
366         fprintf(ftiltl, " %10g", t);
367         if (bKink)
368         {
369             fprintf(fkink, " %10g", t);
370             fprintf(fkinkr, " %10g", t);
371             fprintf(fkinkl, " %10g", t);
372         }
373
374         for (i = 0; i < bun.n; i++)
375         {
376             fprintf(flen, " %6g", bun.len[i]);
377             fprintf(fdist, " %6g", norm(bun.mid[i]));
378             fprintf(fz, " %6g", bun.mid[i][ZZ]);
379             fprintf(ftilt, " %6g", RAD2DEG*acos(bun.dir[i][ZZ]));
380             comp = bun.mid[i][XX]*bun.dir[i][XX]+bun.mid[i][YY]*bun.dir[i][YY];
381             fprintf(ftiltr, " %6g", RAD2DEG*
382                     asin(comp/sqrt(sqr(comp)+sqr(bun.dir[i][ZZ]))));
383             comp = bun.mid[i][YY]*bun.dir[i][XX]-bun.mid[i][XX]*bun.dir[i][YY];
384             fprintf(ftiltl, " %6g", RAD2DEG*
385                     asin(comp/sqrt(sqr(comp)+sqr(bun.dir[i][ZZ]))));
386             if (bKink)
387             {
388                 rvec_sub(bun.end[0][i], bun.end[2][i], va);
389                 rvec_sub(bun.end[2][i], bun.end[1][i], vb);
390                 unitv_no_table(va, va);
391                 unitv_no_table(vb, vb);
392                 fprintf(fkink, " %6g", RAD2DEG*acos(iprod(va, vb)));
393                 cprod(va, vb, vc);
394                 copy_rvec(bun.mid[i], vr);
395                 vr[ZZ] = 0;
396                 unitv_no_table(vr, vr);
397                 fprintf(fkinkr, " %6g", RAD2DEG*asin(iprod(vc, vr)));
398                 vl[XX] = vr[YY];
399                 vl[YY] = -vr[XX];
400                 vl[ZZ] = 0;
401                 fprintf(fkinkl, " %6g", RAD2DEG*asin(iprod(vc, vl)));
402             }
403         }
404         fprintf(flen, "\n");
405         fprintf(fdist, "\n");
406         fprintf(fz, "\n");
407         fprintf(ftilt, "\n");
408         fprintf(ftiltr, "\n");
409         fprintf(ftiltl, "\n");
410         if (bKink)
411         {
412             fprintf(fkink, "\n");
413             fprintf(fkinkr, "\n");
414             fprintf(fkinkl, "\n");
415         }
416         if (fpdb)
417         {
418             dump_axes(fpdb, &fr, &outatoms, &bun);
419         }
420     }
421     while (read_next_frame(oenv, status, &fr));
422     gmx_rmpbc_done(gpbc);
423
424     close_trx(status);
425
426     if (fpdb)
427     {
428         close_trx(fpdb);
429     }
430     gmx_ffclose(flen);
431     gmx_ffclose(fdist);
432     gmx_ffclose(fz);
433     gmx_ffclose(ftilt);
434     gmx_ffclose(ftiltr);
435     gmx_ffclose(ftiltl);
436     if (bKink)
437     {
438         gmx_ffclose(fkink);
439         gmx_ffclose(fkinkr);
440         gmx_ffclose(fkinkl);
441     }
442
443     return 0;
444 }