2 * This file is part of the GROMACS molecular simulation package.
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,2021, 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.
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.
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.
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.
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.
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.
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/utilities.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/pbcutil/rmpbc.h"
53 #include "gromacs/topology/index.h"
54 #include "gromacs/topology/topology.h"
55 #include "gromacs/trajectory/trajectoryframe.h"
56 #include "gromacs/utility/arraysize.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/futil.h"
59 #include "gromacs/utility/smalloc.h"
73 static void rotate_ends(t_bundle* bun, rvec axis, int c0, int c1)
79 for (end = 0; end < bun->nend; end++)
81 for (i = 0; i < bun->n; i++)
83 copy_rvec(bun->end[end][i], tmp);
84 bun->end[end][i][c0] = ax[c1] * tmp[c0] - ax[c0] * tmp[c1];
85 bun->end[end][i][c1] = ax[c0] * tmp[c0] + ax[c1] * tmp[c1];
89 axis[c0] = ax[c1] * tmp[c0] - ax[c0] * tmp[c1];
90 axis[c1] = ax[c0] * tmp[c0] + ax[c1] * tmp[c1];
93 static void calc_axes(rvec x[], t_atom atom[], const int gnx[], int* index[], gmx_bool bRot, t_bundle* bun)
97 rvec axis[MAX_ENDS], cent;
104 for (end = 0; end < bun->nend; end++)
106 for (i = 0; i < bun->n; i++)
108 clear_rvec(bun->end[end][i]);
111 div = gnx[end] / bun->n;
112 for (i = 0; i < gnx[end]; i++)
114 m = atom[index[end][i]].m;
115 for (d = 0; d < DIM; d++)
117 bun->end[end][i / div][d] += m * x[index[end][i]][d];
121 clear_rvec(axis[end]);
122 for (i = 0; i < bun->n; i++)
124 svmul(1.0 / mtot[i], bun->end[end][i], bun->end[end][i]);
125 rvec_inc(axis[end], bun->end[end][i]);
127 svmul(1.0 / bun->n, axis[end], axis[end]);
131 rvec_add(axis[0], axis[1], cent);
132 svmul(0.5, cent, cent);
133 /* center the bundle on the origin */
134 for (end = 0; end < bun->nend; end++)
136 rvec_dec(axis[end], cent);
137 for (i = 0; i < bun->n; i++)
139 rvec_dec(bun->end[end][i], cent);
144 /* rotate the axis parallel to the z-axis */
145 rotate_ends(bun, axis[0], YY, ZZ);
146 rotate_ends(bun, axis[0], XX, ZZ);
148 for (i = 0; i < bun->n; i++)
150 rvec_add(bun->end[0][i], bun->end[1][i], bun->mid[i]);
151 svmul(0.5, bun->mid[i], bun->mid[i]);
152 rvec_sub(bun->end[0][i], bun->end[1][i], bun->dir[i]);
153 bun->len[i] = norm(bun->dir[i]);
154 unitv(bun->dir[i], bun->dir[i]);
158 static void dump_axes(t_trxstatus* status, t_trxframe* fr, t_atoms* outat, t_bundle* bun)
161 static rvec* xout = nullptr;
164 GMX_ASSERT(outat->nr >= bun->n, "");
167 snew(xout, outat->nr);
170 for (i = 0; i < bun->n; i++)
172 copy_rvec(bun->end[0][i], xout[3 * i]);
175 copy_rvec(bun->end[2][i], xout[3 * i + 1]);
179 copy_rvec(bun->mid[i], xout[3 * i + 1]);
181 copy_rvec(bun->end[1][i], xout[3 * i + 2]);
188 frout.natoms = outat->nr;
191 write_trxframe(status, &frout, nullptr);
194 int gmx_bundle(int argc, char* argv[])
196 const char* desc[] = {
197 "[THISMODULE] analyzes bundles of axes. The axes can be for instance",
198 "helix axes. The program reads two index groups and divides both",
199 "of them in [TT]-na[tt] parts. The centers of mass of these parts",
200 "define the tops and bottoms of the axes.",
201 "Several quantities are written to file:",
202 "the axis length, the distance and the z-shift of the axis mid-points",
203 "with respect to the average center of all axes, the total tilt,",
204 "the radial tilt and the lateral tilt with respect to the average axis.",
206 "With options [TT]-ok[tt], [TT]-okr[tt] and [TT]-okl[tt] the total,",
207 "radial and lateral kinks of the axes are plotted. An extra index",
208 "group of kink atoms is required, which is also divided into [TT]-na[tt]",
209 "parts. The kink angle is defined as the angle between the kink-top and",
210 "the bottom-kink vectors.",
212 "With option [TT]-oa[tt] the top, mid (or kink when [TT]-ok[tt] is set)",
213 "and bottom points of each axis",
214 "are written to a [REF].pdb[ref] file each frame. The residue numbers correspond",
215 "to the axis numbers. When viewing this file with Rasmol, use the",
216 "command line option [TT]-nmrpdb[tt], and type [TT]set axis true[tt] to",
217 "display the reference axis."
220 static gmx_bool bZ = FALSE;
221 t_pargs pa[] = { { "-na", FALSE, etINT, { &n }, "Number of axes" },
226 "Use the [IT]z[it]-axis as reference instead of the average axis" } };
227 FILE * flen, *fdist, *fz, *ftilt, *ftiltr, *ftiltl;
228 FILE * fkink = nullptr, *fkinkr = nullptr, *fkinkl = nullptr;
238 char* grpname[MAX_ENDS];
239 /* FIXME: The constness should not be cast away */
240 char * anm = const_cast<char*>("CA"), *rnm = const_cast<char*>("GLY");
241 int i, gnx[MAX_ENDS];
242 int* index[MAX_ENDS];
245 rvec va, vb, vc, vr, vl;
246 gmx_output_env_t* oenv;
247 gmx_rmpbc_t gpbc = nullptr;
249 #define NLEG asize(leg)
251 { efTRX, "-f", nullptr, ffREAD }, { efTPS, nullptr, nullptr, ffREAD },
252 { efNDX, nullptr, nullptr, ffOPTRD }, { efXVG, "-ol", "bun_len", ffWRITE },
253 { efXVG, "-od", "bun_dist", ffWRITE }, { efXVG, "-oz", "bun_z", ffWRITE },
254 { efXVG, "-ot", "bun_tilt", ffWRITE }, { efXVG, "-otr", "bun_tiltr", ffWRITE },
255 { efXVG, "-otl", "bun_tiltl", ffWRITE }, { efXVG, "-ok", "bun_kink", ffOPTWR },
256 { efXVG, "-okr", "bun_kinkr", ffOPTWR }, { efXVG, "-okl", "bun_kinkl", ffOPTWR },
257 { efPDB, "-oa", "axes", ffOPTWR }
259 #define NFILE asize(fnm)
261 if (!parse_common_args(
262 &argc, argv, PCA_CAN_TIME | PCA_TIME_UNIT, NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
267 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &pbcType, &xtop, nullptr, box, TRUE);
269 bKink = opt2bSet("-ok", NFILE, fnm) || opt2bSet("-okr", NFILE, fnm) || opt2bSet("-okl", NFILE, fnm);
279 fprintf(stderr, "Select a group of top and a group of bottom ");
282 fprintf(stderr, "and a group of kink ");
284 fprintf(stderr, "atoms\n");
285 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), bun.nend, gnx, index, grpname);
287 if (n <= 0 || gnx[0] % n || gnx[1] % n || (bKink && gnx[2] % n))
289 gmx_fatal(FARGS, "The size of one of your index groups is not a multiple of n");
303 opt2fn("-ol", NFILE, fnm), "Axis lengths", output_env_get_xvgr_tlabel(oenv), "(nm)", oenv);
304 fdist = xvgropen(opt2fn("-od", NFILE, fnm),
305 "Distance of axis centers",
306 output_env_get_xvgr_tlabel(oenv),
309 fz = xvgropen(opt2fn("-oz", NFILE, fnm),
310 "Z-shift of axis centers",
311 output_env_get_xvgr_tlabel(oenv),
315 opt2fn("-ot", NFILE, fnm), "Axis tilts", output_env_get_xvgr_tlabel(oenv), "(degrees)", oenv);
316 ftiltr = xvgropen(opt2fn("-otr", NFILE, fnm),
318 output_env_get_xvgr_tlabel(oenv),
321 ftiltl = xvgropen(opt2fn("-otl", NFILE, fnm),
322 "Lateral axis tilts",
323 output_env_get_xvgr_tlabel(oenv),
329 fkink = xvgropen(opt2fn("-ok", NFILE, fnm),
331 output_env_get_xvgr_tlabel(oenv),
334 fkinkr = xvgropen(opt2fn("-okr", NFILE, fnm),
335 "Radial kink angles",
336 output_env_get_xvgr_tlabel(oenv),
339 if (output_env_get_print_xvgr_codes(oenv))
341 fprintf(fkinkr, "@ subtitle \"+ = ) ( - = ( )\"\n");
343 fkinkl = xvgropen(opt2fn("-okl", NFILE, fnm),
344 "Lateral kink angles",
345 output_env_get_xvgr_tlabel(oenv),
350 if (opt2bSet("-oa", NFILE, fnm))
352 init_t_atoms(&outatoms, 3 * n, FALSE);
354 for (i = 0; i < 3 * n; i++)
356 outatoms.atomname[i] = &anm;
357 outatoms.atom[i].resind = i / 3;
358 outatoms.resinfo[i / 3].name = &rnm;
359 outatoms.resinfo[i / 3].nr = i / 3 + 1;
360 outatoms.resinfo[i / 3].ic = ' ';
362 fpdb = open_trx(opt2fn("-oa", NFILE, fnm), "w");
369 read_first_frame(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &fr, TRX_NEED_X);
370 gpbc = gmx_rmpbc_init(&top.idef, pbcType, fr.natoms);
374 gmx_rmpbc_trxfr(gpbc, &fr);
375 calc_axes(fr.x, top.atoms.atom, gnx, index, !bZ, &bun);
376 t = output_env_conv_time(oenv, fr.time);
377 fprintf(flen, " %10g", t);
378 fprintf(fdist, " %10g", t);
379 fprintf(fz, " %10g", t);
380 fprintf(ftilt, " %10g", t);
381 fprintf(ftiltr, " %10g", t);
382 fprintf(ftiltl, " %10g", t);
385 fprintf(fkink, " %10g", t);
386 fprintf(fkinkr, " %10g", t);
387 fprintf(fkinkl, " %10g", t);
390 for (i = 0; i < bun.n; i++)
392 fprintf(flen, " %6g", bun.len[i]);
393 fprintf(fdist, " %6g", norm(bun.mid[i]));
394 fprintf(fz, " %6g", bun.mid[i][ZZ]);
395 fprintf(ftilt, " %6g", gmx::c_rad2Deg * acos(bun.dir[i][ZZ]));
396 comp = bun.mid[i][XX] * bun.dir[i][XX] + bun.mid[i][YY] * bun.dir[i][YY];
397 fprintf(ftiltr, " %6g", gmx::c_rad2Deg * std::asin(comp / std::hypot(comp, bun.dir[i][ZZ])));
398 comp = bun.mid[i][YY] * bun.dir[i][XX] - bun.mid[i][XX] * bun.dir[i][YY];
399 fprintf(ftiltl, " %6g", gmx::c_rad2Deg * std::asin(comp / std::hypot(comp, bun.dir[i][ZZ])));
402 rvec_sub(bun.end[0][i], bun.end[2][i], va);
403 rvec_sub(bun.end[2][i], bun.end[1][i], vb);
406 fprintf(fkink, " %6g", gmx::c_rad2Deg * acos(iprod(va, vb)));
408 copy_rvec(bun.mid[i], vr);
411 fprintf(fkinkr, " %6g", gmx::c_rad2Deg * std::asin(iprod(vc, vr)));
415 fprintf(fkinkl, " %6g", gmx::c_rad2Deg * std::asin(iprod(vc, vl)));
419 fprintf(fdist, "\n");
421 fprintf(ftilt, "\n");
422 fprintf(ftiltr, "\n");
423 fprintf(ftiltl, "\n");
426 fprintf(fkink, "\n");
427 fprintf(fkinkr, "\n");
428 fprintf(fkinkl, "\n");
432 dump_axes(fpdb, &fr, &outatoms, &bun);
434 } while (read_next_frame(oenv, status, &fr));
435 gmx_rmpbc_done(gpbc);