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