07564ee2d5b481165ba0445e3ee372c5663984a8
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_bundle.cpp
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,2015,2017,2018 by the GROMACS development team.
7  * Copyright (c) 2019,2020, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #include "gmxpre.h"
39
40 #include <cmath>
41 #include <cstring>
42
43 #include "gromacs/commandline/pargs.h"
44 #include "gromacs/fileio/confio.h"
45 #include "gromacs/fileio/trxio.h"
46 #include "gromacs/fileio/xvgr.h"
47 #include "gromacs/gmxana/gmx_ana.h"
48 #include "gromacs/math/functions.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/topology/topology.h"
54 #include "gromacs/trajectory/trajectoryframe.h"
55 #include "gromacs/utility/arraysize.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/futil.h"
58 #include "gromacs/utility/smalloc.h"
59
60 #define MAX_ENDS 3
61
62 typedef struct
63 {
64     int   n;
65     int   nend;
66     rvec* end[MAX_ENDS];
67     rvec* mid;
68     rvec* dir;
69     real* len;
70 } t_bundle;
71
72 static void rotate_ends(t_bundle* bun, rvec axis, int c0, int c1)
73 {
74     int  end, i;
75     rvec ax, tmp;
76
77     unitv(axis, ax);
78     for (end = 0; end < bun->nend; end++)
79     {
80         for (i = 0; i < bun->n; i++)
81         {
82             copy_rvec(bun->end[end][i], tmp);
83             bun->end[end][i][c0] = ax[c1] * tmp[c0] - ax[c0] * tmp[c1];
84             bun->end[end][i][c1] = ax[c0] * tmp[c0] + ax[c1] * tmp[c1];
85         }
86     }
87     copy_rvec(axis, tmp);
88     axis[c0] = ax[c1] * tmp[c0] - ax[c0] * tmp[c1];
89     axis[c1] = ax[c0] * tmp[c0] + ax[c1] * tmp[c1];
90 }
91
92 static void calc_axes(rvec x[], t_atom atom[], const int gnx[], int* index[], gmx_bool bRot, t_bundle* bun)
93 {
94     int   end, i, div, d;
95     real *mtot, m;
96     rvec  axis[MAX_ENDS], cent;
97
98     snew(mtot, bun->n);
99
100     clear_rvec(axis[0]);
101     clear_rvec(axis[1]);
102
103     for (end = 0; end < bun->nend; end++)
104     {
105         for (i = 0; i < bun->n; i++)
106         {
107             clear_rvec(bun->end[end][i]);
108             mtot[i] = 0;
109         }
110         div = gnx[end] / bun->n;
111         for (i = 0; i < gnx[end]; i++)
112         {
113             m = atom[index[end][i]].m;
114             for (d = 0; d < DIM; d++)
115             {
116                 bun->end[end][i / div][d] += m * x[index[end][i]][d];
117             }
118             mtot[i / div] += m;
119         }
120         clear_rvec(axis[end]);
121         for (i = 0; i < bun->n; i++)
122         {
123             svmul(1.0 / mtot[i], bun->end[end][i], bun->end[end][i]);
124             rvec_inc(axis[end], bun->end[end][i]);
125         }
126         svmul(1.0 / bun->n, axis[end], axis[end]);
127     }
128     sfree(mtot);
129
130     rvec_add(axis[0], axis[1], cent);
131     svmul(0.5, cent, cent);
132     /* center the bundle on the origin */
133     for (end = 0; end < bun->nend; end++)
134     {
135         rvec_dec(axis[end], cent);
136         for (i = 0; i < bun->n; i++)
137         {
138             rvec_dec(bun->end[end][i], cent);
139         }
140     }
141     if (bRot)
142     {
143         /* rotate the axis parallel to the z-axis */
144         rotate_ends(bun, axis[0], YY, ZZ);
145         rotate_ends(bun, axis[0], XX, ZZ);
146     }
147     for (i = 0; i < bun->n; i++)
148     {
149         rvec_add(bun->end[0][i], bun->end[1][i], bun->mid[i]);
150         svmul(0.5, bun->mid[i], bun->mid[i]);
151         rvec_sub(bun->end[0][i], bun->end[1][i], bun->dir[i]);
152         bun->len[i] = norm(bun->dir[i]);
153         unitv(bun->dir[i], bun->dir[i]);
154     }
155 }
156
157 static void dump_axes(t_trxstatus* status, t_trxframe* fr, t_atoms* outat, t_bundle* bun)
158 {
159     t_trxframe   frout;
160     static rvec* xout = nullptr;
161     int          i;
162
163     GMX_ASSERT(outat->nr >= bun->n, "");
164     if (xout == nullptr)
165     {
166         snew(xout, outat->nr);
167     }
168
169     for (i = 0; i < bun->n; i++)
170     {
171         copy_rvec(bun->end[0][i], xout[3 * i]);
172         if (bun->nend >= 3)
173         {
174             copy_rvec(bun->end[2][i], xout[3 * i + 1]);
175         }
176         else
177         {
178             copy_rvec(bun->mid[i], xout[3 * i + 1]);
179         }
180         copy_rvec(bun->end[1][i], xout[3 * i + 2]);
181     }
182     frout        = *fr;
183     frout.bV     = FALSE;
184     frout.bF     = FALSE;
185     frout.bBox   = FALSE;
186     frout.bAtoms = TRUE;
187     frout.natoms = outat->nr;
188     frout.atoms  = outat;
189     frout.x      = xout;
190     write_trxframe(status, &frout, nullptr);
191 }
192
193 int gmx_bundle(int argc, char* argv[])
194 {
195     const char* desc[] = {
196         "[THISMODULE] analyzes bundles of axes. The axes can be for instance",
197         "helix axes. The program reads two index groups and divides both",
198         "of them in [TT]-na[tt] parts. The centers of mass of these parts",
199         "define the tops and bottoms of the axes.",
200         "Several quantities are written to file:",
201         "the axis length, the distance and the z-shift of the axis mid-points",
202         "with respect to the average center of all axes, the total tilt,",
203         "the radial tilt and the lateral tilt with respect to the average axis.",
204         "[PAR]",
205         "With options [TT]-ok[tt], [TT]-okr[tt] and [TT]-okl[tt] the total,",
206         "radial and lateral kinks of the axes are plotted. An extra index",
207         "group of kink atoms is required, which is also divided into [TT]-na[tt]",
208         "parts. The kink angle is defined as the angle between the kink-top and",
209         "the bottom-kink vectors.",
210         "[PAR]",
211         "With option [TT]-oa[tt] the top, mid (or kink when [TT]-ok[tt] is set)",
212         "and bottom points of each axis",
213         "are written to a [REF].pdb[ref] file each frame. The residue numbers correspond",
214         "to the axis numbers. When viewing this file with Rasmol, use the",
215         "command line option [TT]-nmrpdb[tt], and type [TT]set axis true[tt] to",
216         "display the reference axis."
217     };
218     static int      n    = 0;
219     static gmx_bool bZ   = FALSE;
220     t_pargs         pa[] = { { "-na", FALSE, etINT, { &n }, "Number of axes" },
221                      { "-z",
222                        FALSE,
223                        etBOOL,
224                        { &bZ },
225                        "Use the [IT]z[it]-axis as reference instead of the average axis" } };
226     FILE *          flen, *fdist, *fz, *ftilt, *ftiltr, *ftiltl;
227     FILE *          fkink = nullptr, *fkinkr = nullptr, *fkinkl = nullptr;
228     t_trxstatus*    status;
229     t_trxstatus*    fpdb;
230     t_topology      top;
231     PbcType         pbcType;
232     rvec*           xtop;
233     matrix          box;
234     t_trxframe      fr;
235     t_atoms         outatoms;
236     real            t, comp;
237     char*           grpname[MAX_ENDS];
238     /* FIXME: The constness should not be cast away */
239     char *            anm = const_cast<char*>("CA"), *rnm = const_cast<char*>("GLY");
240     int               i, gnx[MAX_ENDS];
241     int*              index[MAX_ENDS];
242     t_bundle          bun;
243     gmx_bool          bKink;
244     rvec              va, vb, vc, vr, vl;
245     gmx_output_env_t* oenv;
246     gmx_rmpbc_t       gpbc = nullptr;
247
248 #define NLEG asize(leg)
249     t_filenm fnm[] = {
250         { efTRX, "-f", nullptr, ffREAD },        { efTPS, nullptr, nullptr, ffREAD },
251         { efNDX, nullptr, nullptr, ffOPTRD },    { efXVG, "-ol", "bun_len", ffWRITE },
252         { efXVG, "-od", "bun_dist", ffWRITE },   { efXVG, "-oz", "bun_z", ffWRITE },
253         { efXVG, "-ot", "bun_tilt", ffWRITE },   { efXVG, "-otr", "bun_tiltr", ffWRITE },
254         { efXVG, "-otl", "bun_tiltl", ffWRITE }, { efXVG, "-ok", "bun_kink", ffOPTWR },
255         { efXVG, "-okr", "bun_kinkr", ffOPTWR }, { efXVG, "-okl", "bun_kinkl", ffOPTWR },
256         { efPDB, "-oa", "axes", ffOPTWR }
257     };
258 #define NFILE asize(fnm)
259
260     if (!parse_common_args(
261                 &argc, argv, PCA_CAN_TIME | PCA_TIME_UNIT, NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
262     {
263         return 0;
264     }
265
266     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &pbcType, &xtop, nullptr, box, TRUE);
267
268     bKink = opt2bSet("-ok", NFILE, fnm) || opt2bSet("-okr", NFILE, fnm) || opt2bSet("-okl", NFILE, fnm);
269     if (bKink)
270     {
271         bun.nend = 3;
272     }
273     else
274     {
275         bun.nend = 2;
276     }
277
278     fprintf(stderr, "Select a group of top and a group of bottom ");
279     if (bKink)
280     {
281         fprintf(stderr, "and a group of kink ");
282     }
283     fprintf(stderr, "atoms\n");
284     get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), bun.nend, gnx, index, grpname);
285
286     if (n <= 0 || gnx[0] % n || gnx[1] % n || (bKink && gnx[2] % n))
287     {
288         gmx_fatal(FARGS, "The size of one of your index groups is not a multiple of n");
289     }
290     bun.n = n;
291     snew(bun.end[0], n);
292     snew(bun.end[1], n);
293     if (bKink)
294     {
295         snew(bun.end[2], n);
296     }
297     snew(bun.mid, n);
298     snew(bun.dir, n);
299     snew(bun.len, n);
300
301     flen = xvgropen(
302             opt2fn("-ol", NFILE, fnm), "Axis lengths", output_env_get_xvgr_tlabel(oenv), "(nm)", oenv);
303     fdist = xvgropen(opt2fn("-od", NFILE, fnm),
304                      "Distance of axis centers",
305                      output_env_get_xvgr_tlabel(oenv),
306                      "(nm)",
307                      oenv);
308     fz    = xvgropen(opt2fn("-oz", NFILE, fnm),
309                   "Z-shift of axis centers",
310                   output_env_get_xvgr_tlabel(oenv),
311                   "(nm)",
312                   oenv);
313     ftilt = xvgropen(
314             opt2fn("-ot", NFILE, fnm), "Axis tilts", output_env_get_xvgr_tlabel(oenv), "(degrees)", oenv);
315     ftiltr = xvgropen(opt2fn("-otr", NFILE, fnm),
316                       "Radial axis tilts",
317                       output_env_get_xvgr_tlabel(oenv),
318                       "(degrees)",
319                       oenv);
320     ftiltl = xvgropen(opt2fn("-otl", NFILE, fnm),
321                       "Lateral axis tilts",
322                       output_env_get_xvgr_tlabel(oenv),
323                       "(degrees)",
324                       oenv);
325
326     if (bKink)
327     {
328         fkink  = xvgropen(opt2fn("-ok", NFILE, fnm),
329                          "Kink angles",
330                          output_env_get_xvgr_tlabel(oenv),
331                          "(degrees)",
332                          oenv);
333         fkinkr = xvgropen(opt2fn("-okr", NFILE, fnm),
334                           "Radial kink angles",
335                           output_env_get_xvgr_tlabel(oenv),
336                           "(degrees)",
337                           oenv);
338         if (output_env_get_print_xvgr_codes(oenv))
339         {
340             fprintf(fkinkr, "@ subtitle \"+ = ) (   - = ( )\"\n");
341         }
342         fkinkl = xvgropen(opt2fn("-okl", NFILE, fnm),
343                           "Lateral kink angles",
344                           output_env_get_xvgr_tlabel(oenv),
345                           "(degrees)",
346                           oenv);
347     }
348
349     if (opt2bSet("-oa", NFILE, fnm))
350     {
351         init_t_atoms(&outatoms, 3 * n, FALSE);
352         outatoms.nr = 3 * n;
353         for (i = 0; i < 3 * n; i++)
354         {
355             outatoms.atomname[i]         = &anm;
356             outatoms.atom[i].resind      = i / 3;
357             outatoms.resinfo[i / 3].name = &rnm;
358             outatoms.resinfo[i / 3].nr   = i / 3 + 1;
359             outatoms.resinfo[i / 3].ic   = ' ';
360         }
361         fpdb = open_trx(opt2fn("-oa", NFILE, fnm), "w");
362     }
363     else
364     {
365         fpdb = nullptr;
366     }
367
368     read_first_frame(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &fr, TRX_NEED_X);
369     gpbc = gmx_rmpbc_init(&top.idef, pbcType, fr.natoms);
370
371     do
372     {
373         gmx_rmpbc_trxfr(gpbc, &fr);
374         calc_axes(fr.x, top.atoms.atom, gnx, index, !bZ, &bun);
375         t = output_env_conv_time(oenv, fr.time);
376         fprintf(flen, " %10g", t);
377         fprintf(fdist, " %10g", t);
378         fprintf(fz, " %10g", t);
379         fprintf(ftilt, " %10g", t);
380         fprintf(ftiltr, " %10g", t);
381         fprintf(ftiltl, " %10g", t);
382         if (bKink)
383         {
384             fprintf(fkink, " %10g", t);
385             fprintf(fkinkr, " %10g", t);
386             fprintf(fkinkl, " %10g", t);
387         }
388
389         for (i = 0; i < bun.n; i++)
390         {
391             fprintf(flen, " %6g", bun.len[i]);
392             fprintf(fdist, " %6g", norm(bun.mid[i]));
393             fprintf(fz, " %6g", bun.mid[i][ZZ]);
394             fprintf(ftilt, " %6g", RAD2DEG * acos(bun.dir[i][ZZ]));
395             comp = bun.mid[i][XX] * bun.dir[i][XX] + bun.mid[i][YY] * bun.dir[i][YY];
396             fprintf(ftiltr, " %6g", RAD2DEG * std::asin(comp / std::hypot(comp, bun.dir[i][ZZ])));
397             comp = bun.mid[i][YY] * bun.dir[i][XX] - bun.mid[i][XX] * bun.dir[i][YY];
398             fprintf(ftiltl, " %6g", RAD2DEG * std::asin(comp / std::hypot(comp, bun.dir[i][ZZ])));
399             if (bKink)
400             {
401                 rvec_sub(bun.end[0][i], bun.end[2][i], va);
402                 rvec_sub(bun.end[2][i], bun.end[1][i], vb);
403                 unitv(va, va);
404                 unitv(vb, vb);
405                 fprintf(fkink, " %6g", RAD2DEG * acos(iprod(va, vb)));
406                 cprod(va, vb, vc);
407                 copy_rvec(bun.mid[i], vr);
408                 vr[ZZ] = 0;
409                 unitv(vr, vr);
410                 fprintf(fkinkr, " %6g", RAD2DEG * std::asin(iprod(vc, vr)));
411                 vl[XX] = vr[YY];
412                 vl[YY] = -vr[XX];
413                 vl[ZZ] = 0;
414                 fprintf(fkinkl, " %6g", RAD2DEG * std::asin(iprod(vc, vl)));
415             }
416         }
417         fprintf(flen, "\n");
418         fprintf(fdist, "\n");
419         fprintf(fz, "\n");
420         fprintf(ftilt, "\n");
421         fprintf(ftiltr, "\n");
422         fprintf(ftiltl, "\n");
423         if (bKink)
424         {
425             fprintf(fkink, "\n");
426             fprintf(fkinkr, "\n");
427             fprintf(fkinkl, "\n");
428         }
429         if (fpdb)
430         {
431             dump_axes(fpdb, &fr, &outatoms, &bun);
432         }
433     } while (read_next_frame(oenv, status, &fr));
434     gmx_rmpbc_done(gpbc);
435
436     close_trx(status);
437
438     if (fpdb)
439     {
440         close_trx(fpdb);
441     }
442     xvgrclose(flen);
443     xvgrclose(fdist);
444     xvgrclose(fz);
445     xvgrclose(ftilt);
446     xvgrclose(ftiltr);
447     xvgrclose(ftiltl);
448     if (bKink)
449     {
450         xvgrclose(fkink);
451         xvgrclose(fkinkr);
452         xvgrclose(fkinkl);
453     }
454
455     return 0;
456 }