3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Green Red Orange Magenta Azure Cyan Skyblue
69 static void rotate_ends(t_bundle *bun, rvec axis, int c0, int c1)
75 for (end = 0; end < bun->nend; end++)
77 for (i = 0; i < bun->n; i++)
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];
85 axis[c0] = ax[c1]*tmp[c0] - ax[c0]*tmp[c1];
86 axis[c1] = ax[c0]*tmp[c0] + ax[c1]*tmp[c1];
89 static void calc_axes(rvec x[], t_atom atom[], int gnx[], atom_id *index[],
90 gmx_bool bRot, t_bundle *bun)
94 rvec axis[MAX_ENDS], cent;
98 for (end = 0; end < bun->nend; end++)
100 for (i = 0; i < bun->n; i++)
102 clear_rvec(bun->end[end][i]);
105 div = gnx[end]/bun->n;
106 for (i = 0; i < gnx[end]; i++)
108 m = atom[index[end][i]].m;
109 for (d = 0; d < DIM; d++)
111 bun->end[end][i/div][d] += m*x[index[end][i]][d];
115 clear_rvec(axis[end]);
116 for (i = 0; i < bun->n; i++)
118 svmul(1.0/mtot[i], bun->end[end][i], bun->end[end][i]);
119 rvec_inc(axis[end], bun->end[end][i]);
121 svmul(1.0/bun->n, axis[end], axis[end]);
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++)
130 rvec_dec(axis[end], cent);
131 for (i = 0; i < bun->n; i++)
133 rvec_dec(bun->end[end][i], cent);
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);
142 for (i = 0; i < bun->n; i++)
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]);
152 static void dump_axes(t_trxstatus *status, t_trxframe *fr, t_atoms *outat,
156 static rvec *xout = NULL;
161 snew(xout, outat->nr);
164 for (i = 0; i < bun->n; i++)
166 copy_rvec(bun->end[0][i], xout[3*i]);
169 copy_rvec(bun->end[2][i], xout[3*i+1]);
173 copy_rvec(bun->mid[i], xout[3*i+1]);
175 copy_rvec(bun->end[1][i], xout[3*i+2]);
182 frout.natoms = outat->nr;
185 write_trxframe(status, &frout, NULL);
188 int gmx_bundle(int argc, char *argv[])
190 const char *desc[] = {
191 "[TT]g_bundle[tt] 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.",
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.",
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."
214 static gmx_bool bZ = FALSE;
216 { "-na", FALSE, etINT, {&n},
218 { "-z", FALSE, etBOOL, {&bZ},
219 "Use the [IT]z[it]-axis as reference instead of the average axis" }
221 FILE *out, *flen, *fdist, *fz, *ftilt, *ftiltr, *ftiltl;
222 FILE *fkink = NULL, *fkinkr = NULL, *fkinkl = NULL;
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];
240 rvec va, vb, vc, vr, vl;
242 gmx_rmpbc_t gpbc = NULL;
244 #define NLEG asize(leg)
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 }
260 #define NFILE asize(fnm)
262 parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_TIME_UNIT | PCA_BE_NICE,
263 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv);
265 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xtop, NULL, box, TRUE);
267 bKink = opt2bSet("-ok", NFILE, fnm) || opt2bSet("-okr", NFILE, fnm)
268 || opt2bSet("-okl", NFILE, fnm);
278 fprintf(stderr, "Select a group of top and a group of bottom ");
281 fprintf(stderr, "and a group of kink ");
283 fprintf(stderr, "atoms\n");
284 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), bun.nend,
285 gnx, index, grpname);
287 if (n <= 0 || gnx[0] % n || gnx[1] % n || (bKink && gnx[2] % n))
290 "The size of one of your index groups is not a multiple of n");
303 flen = xvgropen(opt2fn("-ol", NFILE, fnm), "Axis lengths",
304 output_env_get_xvgr_tlabel(oenv), "(nm)", oenv);
305 fdist = xvgropen(opt2fn("-od", NFILE, fnm), "Distance of axis centers",
306 output_env_get_xvgr_tlabel(oenv), "(nm)", oenv);
307 fz = xvgropen(opt2fn("-oz", NFILE, fnm), "Z-shift of axis centers",
308 output_env_get_xvgr_tlabel(oenv), "(nm)", oenv);
309 ftilt = xvgropen(opt2fn("-ot", NFILE, fnm), "Axis tilts",
310 output_env_get_xvgr_tlabel(oenv), "(degrees)", oenv);
311 ftiltr = xvgropen(opt2fn("-otr", NFILE, fnm), "Radial axis tilts",
312 output_env_get_xvgr_tlabel(oenv), "(degrees)", oenv);
313 ftiltl = xvgropen(opt2fn("-otl", NFILE, fnm), "Lateral axis tilts",
314 output_env_get_xvgr_tlabel(oenv), "(degrees)", oenv);
318 fkink = xvgropen(opt2fn("-ok", NFILE, fnm), "Kink angles",
319 output_env_get_xvgr_tlabel(oenv), "(degrees)", oenv);
320 fkinkr = xvgropen(opt2fn("-okr", NFILE, fnm), "Radial kink angles",
321 output_env_get_xvgr_tlabel(oenv), "(degrees)", oenv);
322 if (output_env_get_print_xvgr_codes(oenv))
324 fprintf(fkinkr, "@ subtitle \"+ = ) ( - = ( )\"\n");
326 fkinkl = xvgropen(opt2fn("-okl", NFILE, fnm), "Lateral kink angles",
327 output_env_get_xvgr_tlabel(oenv), "(degrees)", oenv);
330 if (opt2bSet("-oa", NFILE, fnm))
332 init_t_atoms(&outatoms, 3*n, FALSE);
334 for (i = 0; i < 3*n; i++)
336 outatoms.atomname[i] = &anm;
337 outatoms.atom[i].resind = i/3;
338 outatoms.resinfo[i/3].name = &rnm;
339 outatoms.resinfo[i/3].nr = i/3 + 1;
340 outatoms.resinfo[i/3].ic = ' ';
342 fpdb = open_trx(opt2fn("-oa", NFILE, fnm), "w");
349 read_first_frame(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &fr, TRX_NEED_X);
350 gpbc = gmx_rmpbc_init(&top.idef, ePBC, fr.natoms, fr.box);
354 gmx_rmpbc_trxfr(gpbc, &fr);
355 calc_axes(fr.x, top.atoms.atom, gnx, index, !bZ, &bun);
356 t = output_env_conv_time(oenv, fr.time);
357 fprintf(flen, " %10g", t);
358 fprintf(fdist, " %10g", t);
359 fprintf(fz, " %10g", t);
360 fprintf(ftilt, " %10g", t);
361 fprintf(ftiltr, " %10g", t);
362 fprintf(ftiltl, " %10g", t);
365 fprintf(fkink, " %10g", t);
366 fprintf(fkinkr, " %10g", t);
367 fprintf(fkinkl, " %10g", t);
370 for (i = 0; i < bun.n; i++)
372 fprintf(flen, " %6g", bun.len[i]);
373 fprintf(fdist, " %6g", norm(bun.mid[i]));
374 fprintf(fz, " %6g", bun.mid[i][ZZ]);
375 fprintf(ftilt, " %6g", RAD2DEG*acos(bun.dir[i][ZZ]));
376 comp = bun.mid[i][XX]*bun.dir[i][XX]+bun.mid[i][YY]*bun.dir[i][YY];
377 fprintf(ftiltr, " %6g", RAD2DEG*
378 asin(comp/sqrt(sqr(comp)+sqr(bun.dir[i][ZZ]))));
379 comp = bun.mid[i][YY]*bun.dir[i][XX]-bun.mid[i][XX]*bun.dir[i][YY];
380 fprintf(ftiltl, " %6g", RAD2DEG*
381 asin(comp/sqrt(sqr(comp)+sqr(bun.dir[i][ZZ]))));
384 rvec_sub(bun.end[0][i], bun.end[2][i], va);
385 rvec_sub(bun.end[2][i], bun.end[1][i], vb);
386 unitv_no_table(va, va);
387 unitv_no_table(vb, vb);
388 fprintf(fkink, " %6g", RAD2DEG*acos(iprod(va, vb)));
390 copy_rvec(bun.mid[i], vr);
392 unitv_no_table(vr, vr);
393 fprintf(fkinkr, " %6g", RAD2DEG*asin(iprod(vc, vr)));
397 fprintf(fkinkl, " %6g", RAD2DEG*asin(iprod(vc, vl)));
401 fprintf(fdist, "\n");
403 fprintf(ftilt, "\n");
404 fprintf(ftiltr, "\n");
405 fprintf(ftiltl, "\n");
408 fprintf(fkink, "\n");
409 fprintf(fkinkr, "\n");
410 fprintf(fkinkl, "\n");
414 dump_axes(fpdb, &fr, &outatoms, &bun);
417 while (read_next_frame(oenv, status, &fr));
418 gmx_rmpbc_done(gpbc);