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
44 #include "gmx_fatal.h"
52 static void rot_conf(t_atoms *atoms, rvec x[], rvec v[], real trans, real angle,
53 rvec head, rvec tail, int isize, atom_id index[],
54 rvec xout[], rvec vout[])
57 real theta, phi, arrow_len;
58 mat4 Rx, Ry, Rz, Rinvy, Rinvz, Mtot;
59 mat4 temp1, temp2, temp3;
63 rvec_sub(tail, head, arrow);
64 arrow_len = norm(arrow);
67 fprintf(debug, "Arrow vector: %10.4f %10.4f %10.4f\n",
68 arrow[XX], arrow[YY], arrow[ZZ]);
69 fprintf(debug, "Effective translation %g nm\n", trans);
73 gmx_fatal(FARGS, "Arrow vector not given");
76 /* Copy all aoms to output */
77 for (i = 0; (i < atoms->nr); i++)
79 copy_rvec(x[i], xout[i]);
80 copy_rvec(v[i], vout[i]);
83 /* Compute center of mass and move atoms there */
85 for (i = 0; (i < isize); i++)
87 rvec_inc(xcm, x[index[i]]);
89 for (i = 0; (i < DIM); i++)
95 fprintf(debug, "Center of mass: %10.4f %10.4f %10.4f\n",
96 xcm[XX], xcm[YY], xcm[ZZ]);
98 for (i = 0; (i < isize); i++)
100 rvec_sub(x[index[i]], xcm, xout[index[i]]);
103 /* Compute theta and phi that describe the arrow */
104 theta = acos(arrow[ZZ]/arrow_len);
105 phi = atan2(arrow[YY]/arrow_len, arrow[XX]/arrow_len);
108 fprintf(debug, "Phi = %.1f, Theta = %.1f\n", RAD2DEG*phi, RAD2DEG*theta);
111 /* Now the total rotation matrix: */
112 /* Rotate a couple of times */
113 rotate(ZZ, -phi, Rz);
114 rotate(YY, M_PI/2-theta, Ry);
115 rotate(XX, angle*DEG2RAD, Rx);
117 rotate(YY, theta-M_PI/2, Rinvy);
118 rotate(ZZ, phi, Rinvz);
120 mult_matrix(temp1, Ry, Rz);
121 mult_matrix(temp2, Rinvy, Rx);
122 mult_matrix(temp3, temp2, temp1);
123 mult_matrix(Mtot, Rinvz, temp3);
125 print_m4(debug, "Rz", Rz);
126 print_m4(debug, "Ry", Ry);
127 print_m4(debug, "Rx", Rx);
128 print_m4(debug, "Rinvy", Rinvy);
129 print_m4(debug, "Rinvz", Rinvz);
130 print_m4(debug, "Mtot", Mtot);
132 for (i = 0; (i < isize); i++)
135 m4_op(Mtot, xout[ai], xv);
136 rvec_add(xv, xcm, xout[ai]);
137 m4_op(Mtot, v[ai], xv);
138 copy_rvec(xv, vout[ai]);
142 int gmx_dyndom(int argc, char *argv[])
144 const char *desc[] = {
145 "[TT]g_dyndom[tt] reads a [TT].pdb[tt] file output from DynDom",
146 "(http://www.cmp.uea.ac.uk/dyndom/).",
147 "It reads the coordinates, the coordinates of the rotation axis,",
148 "and an index file containing the domains.",
149 "Furthermore, it takes the first and last atom of the arrow file",
150 "as command line arguments (head and tail) and",
151 "finally it takes the translation vector (given in DynDom info file)",
152 "and the angle of rotation (also as command line arguments). If the angle",
153 "determined by DynDom is given, one should be able to recover the",
154 "second structure used for generating the DynDom output.",
155 "Because of limited numerical accuracy this should be verified by",
156 "computing an all-atom RMSD (using [TT]g_confrms[tt]) rather than by file",
157 "comparison (using diff).[PAR]",
158 "The purpose of this program is to interpolate and extrapolate the",
159 "rotation as found by DynDom. As a result unphysical structures with",
160 "long or short bonds, or overlapping atoms may be produced. Visual",
161 "inspection, and energy minimization may be necessary to",
162 "validate the structure."
164 static real trans0 = 0;
165 static rvec head = { 0, 0, 0 };
166 static rvec tail = { 0, 0, 0 };
167 static real angle0 = 0, angle1 = 0, maxangle = 0;
168 static int label = 0, nframes = 11;
170 { "-firstangle", FALSE, etREAL, {&angle0},
171 "Angle of rotation about rotation vector" },
172 { "-lastangle", FALSE, etREAL, {&angle1},
173 "Angle of rotation about rotation vector" },
174 { "-nframe", FALSE, etINT, {&nframes},
175 "Number of steps on the pathway" },
176 { "-maxangle", FALSE, etREAL, {&maxangle},
177 "DymDom dtermined angle of rotation about rotation vector" },
178 { "-trans", FALSE, etREAL, {&trans0},
179 "Translation (Angstrom) along rotation vector (see DynDom info file)" },
180 { "-head", FALSE, etRVEC, {head},
181 "First atom of the arrow vector" },
182 { "-tail", FALSE, etRVEC, {tail},
183 "Last atom of the arrow vector" }
185 int i, j, natoms, isize;
187 atom_id *index = NULL, *index_all;
188 char title[256], *grpname;
191 rvec *x, *v, *xout, *vout;
196 { efPDB, "-f", "dyndom", ffREAD },
197 { efTRO, "-o", "rotated", ffWRITE },
198 { efNDX, "-n", "domains", ffREAD }
200 #define NFILE asize(fnm)
202 if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
203 asize(desc), desc, 0, NULL, &oenv))
208 get_stx_coordnum (opt2fn("-f", NFILE, fnm), &natoms);
209 init_t_atoms(&atoms, natoms, TRUE);
212 read_stx_conf(opt2fn("-f", NFILE, fnm), title, &atoms, x, v, NULL, box);
216 printf("Select group to rotate:\n");
217 rd_index(ftp2fn(efNDX, NFILE, fnm), 1, &isize, &index, &grpname);
218 printf("Going to rotate %s containg %d atoms\n", grpname, isize);
220 snew(index_all, atoms.nr);
221 for (i = 0; (i < atoms.nr); i++)
226 status = open_trx(opt2fn("-o", NFILE, fnm), "w");
229 for (i = 0; (i < nframes); i++, label++)
231 angle = angle0 + (i*(angle1-angle0))/(nframes-1);
232 trans = trans0*0.1*angle/maxangle;
233 printf("Frame: %2d (label %c), angle: %8.3f deg., trans: %8.3f nm\n",
234 i, label, angle, trans);
235 rot_conf(&atoms, x, v, trans, angle, head, tail, isize, index, xout, vout);
241 for (j = 0; (j < atoms.nr); j++)
243 atoms.resinfo[atoms.atom[j].resind].chainid = label;
246 write_trx(status, atoms.nr, index_all, &atoms, i, angle, box, xout, vout, NULL);