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