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