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, 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.
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.
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.
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.
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.
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.
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"
70 static void rotate_ends(t_bundle
*bun
, rvec axis
, int c0
, int c1
)
76 for (end
= 0; end
< bun
->nend
; end
++)
78 for (i
= 0; i
< bun
->n
; i
++)
80 copy_rvec(bun
->end
[end
][i
], tmp
);
81 bun
->end
[end
][i
][c0
] = ax
[c1
]*tmp
[c0
] - ax
[c0
]*tmp
[c1
];
82 bun
->end
[end
][i
][c1
] = ax
[c0
]*tmp
[c0
] + ax
[c1
]*tmp
[c1
];
86 axis
[c0
] = ax
[c1
]*tmp
[c0
] - ax
[c0
]*tmp
[c1
];
87 axis
[c1
] = ax
[c0
]*tmp
[c0
] + ax
[c1
]*tmp
[c1
];
90 static void calc_axes(rvec x
[], t_atom atom
[], const int gnx
[], int *index
[],
91 gmx_bool bRot
, t_bundle
*bun
)
95 rvec axis
[MAX_ENDS
], cent
;
102 for (end
= 0; end
< bun
->nend
; end
++)
104 for (i
= 0; i
< bun
->n
; i
++)
106 clear_rvec(bun
->end
[end
][i
]);
109 div
= gnx
[end
]/bun
->n
;
110 for (i
= 0; i
< gnx
[end
]; i
++)
112 m
= atom
[index
[end
][i
]].m
;
113 for (d
= 0; d
< DIM
; d
++)
115 bun
->end
[end
][i
/div
][d
] += m
*x
[index
[end
][i
]][d
];
119 clear_rvec(axis
[end
]);
120 for (i
= 0; i
< bun
->n
; i
++)
122 svmul(1.0/mtot
[i
], bun
->end
[end
][i
], bun
->end
[end
][i
]);
123 rvec_inc(axis
[end
], bun
->end
[end
][i
]);
125 svmul(1.0/bun
->n
, axis
[end
], axis
[end
]);
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
++)
134 rvec_dec(axis
[end
], cent
);
135 for (i
= 0; i
< bun
->n
; i
++)
137 rvec_dec(bun
->end
[end
][i
], cent
);
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
);
146 for (i
= 0; i
< bun
->n
; i
++)
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
]);
156 static void dump_axes(t_trxstatus
*status
, t_trxframe
*fr
, t_atoms
*outat
,
160 static rvec
*xout
= nullptr;
163 GMX_ASSERT(outat
->nr
>= bun
->n
, "");
166 snew(xout
, outat
->nr
);
169 for (i
= 0; i
< bun
->n
; i
++)
171 copy_rvec(bun
->end
[0][i
], xout
[3*i
]);
174 copy_rvec(bun
->end
[2][i
], xout
[3*i
+1]);
178 copy_rvec(bun
->mid
[i
], xout
[3*i
+1]);
180 copy_rvec(bun
->end
[1][i
], xout
[3*i
+2]);
187 frout
.natoms
= outat
->nr
;
190 write_trxframe(status
, &frout
, nullptr);
193 int gmx_bundle(int argc
, char *argv
[])
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.",
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.",
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."
219 static gmx_bool bZ
= FALSE
;
221 { "-na", FALSE
, etINT
, {&n
},
223 { "-z", FALSE
, etBOOL
, {&bZ
},
224 "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;
237 char *grpname
[MAX_ENDS
];
238 /* FIXME: The constness should not be cast away */
239 char *anm
= (char *)"CA", *rnm
= (char *)"GLY";
240 int i
, gnx
[MAX_ENDS
];
241 int *index
[MAX_ENDS
];
244 rvec va
, vb
, vc
, vr
, vl
;
245 gmx_output_env_t
*oenv
;
246 gmx_rmpbc_t gpbc
= nullptr;
248 #define NLEG asize(leg)
250 { efTRX
, "-f", nullptr, ffREAD
},
251 { efTPS
, nullptr, nullptr, ffREAD
},
252 { efNDX
, nullptr, nullptr, ffOPTRD
},
253 { efXVG
, "-ol", "bun_len", ffWRITE
},
254 { efXVG
, "-od", "bun_dist", ffWRITE
},
255 { efXVG
, "-oz", "bun_z", ffWRITE
},
256 { efXVG
, "-ot", "bun_tilt", ffWRITE
},
257 { efXVG
, "-otr", "bun_tiltr", ffWRITE
},
258 { efXVG
, "-otl", "bun_tiltl", ffWRITE
},
259 { efXVG
, "-ok", "bun_kink", ffOPTWR
},
260 { efXVG
, "-okr", "bun_kinkr", ffOPTWR
},
261 { efXVG
, "-okl", "bun_kinkl", ffOPTWR
},
262 { efPDB
, "-oa", "axes", ffOPTWR
}
264 #define NFILE asize(fnm)
266 if (!parse_common_args(&argc
, argv
, PCA_CAN_TIME
| PCA_TIME_UNIT
,
267 NFILE
, fnm
, asize(pa
), pa
, asize(desc
), desc
, 0, nullptr, &oenv
))
272 read_tps_conf(ftp2fn(efTPS
, NFILE
, fnm
), &top
, &ePBC
, &xtop
, nullptr, box
, TRUE
);
274 bKink
= opt2bSet("-ok", NFILE
, fnm
) || opt2bSet("-okr", NFILE
, fnm
)
275 || opt2bSet("-okl", NFILE
, fnm
);
285 fprintf(stderr
, "Select a group of top and a group of bottom ");
288 fprintf(stderr
, "and a group of kink ");
290 fprintf(stderr
, "atoms\n");
291 get_index(&top
.atoms
, ftp2fn_null(efNDX
, NFILE
, fnm
), bun
.nend
,
292 gnx
, index
, grpname
);
294 if (n
<= 0 || gnx
[0] % n
|| gnx
[1] % n
|| (bKink
&& gnx
[2] % n
))
297 "The size of one of your index groups is not a multiple of n");
310 flen
= xvgropen(opt2fn("-ol", NFILE
, fnm
), "Axis lengths",
311 output_env_get_xvgr_tlabel(oenv
), "(nm)", oenv
);
312 fdist
= xvgropen(opt2fn("-od", NFILE
, fnm
), "Distance of axis centers",
313 output_env_get_xvgr_tlabel(oenv
), "(nm)", oenv
);
314 fz
= xvgropen(opt2fn("-oz", NFILE
, fnm
), "Z-shift of axis centers",
315 output_env_get_xvgr_tlabel(oenv
), "(nm)", oenv
);
316 ftilt
= xvgropen(opt2fn("-ot", NFILE
, fnm
), "Axis tilts",
317 output_env_get_xvgr_tlabel(oenv
), "(degrees)", oenv
);
318 ftiltr
= xvgropen(opt2fn("-otr", NFILE
, fnm
), "Radial axis tilts",
319 output_env_get_xvgr_tlabel(oenv
), "(degrees)", oenv
);
320 ftiltl
= xvgropen(opt2fn("-otl", NFILE
, fnm
), "Lateral axis tilts",
321 output_env_get_xvgr_tlabel(oenv
), "(degrees)", oenv
);
325 fkink
= xvgropen(opt2fn("-ok", NFILE
, fnm
), "Kink angles",
326 output_env_get_xvgr_tlabel(oenv
), "(degrees)", oenv
);
327 fkinkr
= xvgropen(opt2fn("-okr", NFILE
, fnm
), "Radial kink angles",
328 output_env_get_xvgr_tlabel(oenv
), "(degrees)", oenv
);
329 if (output_env_get_print_xvgr_codes(oenv
))
331 fprintf(fkinkr
, "@ subtitle \"+ = ) ( - = ( )\"\n");
333 fkinkl
= xvgropen(opt2fn("-okl", NFILE
, fnm
), "Lateral kink angles",
334 output_env_get_xvgr_tlabel(oenv
), "(degrees)", oenv
);
337 if (opt2bSet("-oa", NFILE
, fnm
))
339 init_t_atoms(&outatoms
, 3*n
, FALSE
);
341 for (i
= 0; i
< 3*n
; i
++)
343 outatoms
.atomname
[i
] = &anm
;
344 outatoms
.atom
[i
].resind
= i
/3;
345 outatoms
.resinfo
[i
/3].name
= &rnm
;
346 outatoms
.resinfo
[i
/3].nr
= i
/3 + 1;
347 outatoms
.resinfo
[i
/3].ic
= ' ';
349 fpdb
= open_trx(opt2fn("-oa", NFILE
, fnm
), "w");
356 read_first_frame(oenv
, &status
, ftp2fn(efTRX
, NFILE
, fnm
), &fr
, TRX_NEED_X
);
357 gpbc
= gmx_rmpbc_init(&top
.idef
, ePBC
, fr
.natoms
);
361 gmx_rmpbc_trxfr(gpbc
, &fr
);
362 calc_axes(fr
.x
, top
.atoms
.atom
, gnx
, index
, !bZ
, &bun
);
363 t
= output_env_conv_time(oenv
, fr
.time
);
364 fprintf(flen
, " %10g", t
);
365 fprintf(fdist
, " %10g", t
);
366 fprintf(fz
, " %10g", t
);
367 fprintf(ftilt
, " %10g", t
);
368 fprintf(ftiltr
, " %10g", t
);
369 fprintf(ftiltl
, " %10g", t
);
372 fprintf(fkink
, " %10g", t
);
373 fprintf(fkinkr
, " %10g", t
);
374 fprintf(fkinkl
, " %10g", t
);
377 for (i
= 0; i
< bun
.n
; i
++)
379 fprintf(flen
, " %6g", bun
.len
[i
]);
380 fprintf(fdist
, " %6g", norm(bun
.mid
[i
]));
381 fprintf(fz
, " %6g", bun
.mid
[i
][ZZ
]);
382 fprintf(ftilt
, " %6g", RAD2DEG
*acos(bun
.dir
[i
][ZZ
]));
383 comp
= bun
.mid
[i
][XX
]*bun
.dir
[i
][XX
]+bun
.mid
[i
][YY
]*bun
.dir
[i
][YY
];
384 fprintf(ftiltr
, " %6g", RAD2DEG
*
385 std::asin(comp
/std::hypot(comp
, bun
.dir
[i
][ZZ
])));
386 comp
= bun
.mid
[i
][YY
]*bun
.dir
[i
][XX
]-bun
.mid
[i
][XX
]*bun
.dir
[i
][YY
];
387 fprintf(ftiltl
, " %6g", RAD2DEG
*
388 std::asin(comp
/std::hypot(comp
, bun
.dir
[i
][ZZ
])));
391 rvec_sub(bun
.end
[0][i
], bun
.end
[2][i
], va
);
392 rvec_sub(bun
.end
[2][i
], bun
.end
[1][i
], vb
);
395 fprintf(fkink
, " %6g", RAD2DEG
*acos(iprod(va
, vb
)));
397 copy_rvec(bun
.mid
[i
], vr
);
400 fprintf(fkinkr
, " %6g", RAD2DEG
*std::asin(iprod(vc
, vr
)));
404 fprintf(fkinkl
, " %6g", RAD2DEG
*std::asin(iprod(vc
, vl
)));
408 fprintf(fdist
, "\n");
410 fprintf(ftilt
, "\n");
411 fprintf(ftiltr
, "\n");
412 fprintf(ftiltl
, "\n");
415 fprintf(fkink
, "\n");
416 fprintf(fkinkr
, "\n");
417 fprintf(fkinkl
, "\n");
421 dump_axes(fpdb
, &fr
, &outatoms
, &bun
);
424 while (read_next_frame(oenv
, status
, &fr
));
425 gmx_rmpbc_done(gpbc
);