1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * Green Red Orange Magenta Azure Cyan Skyblue
62 static void get_refx(output_env_t oenv, const char *trxfn, int nfitdim, int skip,
64 gmx_bool bMW, t_topology *top, int ePBC, rvec *x_ref)
66 int natoms, nfr_all, nfr, i, j, a, r, c, min_fr;
69 double tot_mass, msd, *srmsd, min_srmsd, srmsd_tot;
74 gmx_rmpbc_t gpbc = NULL;
81 natoms = read_first_x(oenv, &status, trxfn, &ti[nfr], &x, box);
85 for (a = 0; a < gnx; a++)
87 if (index[a] >= natoms)
89 gmx_fatal(FARGS, "Atom index (%d) is larger than the number of atoms in the trajecory (%d)", index[a]+1, natoms);
91 w_rls[a] = (bMW ? top->atoms.atom[index[a]].m : 1.0);
94 gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms, box);
98 if (nfr_all % skip == 0)
100 gmx_rmpbc(gpbc, natoms, box, x);
102 for (i = 0; i < gnx; i++)
104 copy_rvec(x[index[i]], xi[nfr][i]);
106 reset_x(gnx, NULL, gnx, NULL, xi[nfr], w_rls);
116 while (read_next_x(oenv, status, &ti[nfr], natoms, x, box));
120 gmx_rmpbc_done(gpbc);
123 for (i = 0; i < nfr; i++)
125 printf("\rProcessing frame %d of %d", i, nfr);
126 for (j = i+1; j < nfr; j++)
128 calc_fit_R(nfitdim, gnx, w_rls, xi[i], xi[j], R);
131 for (a = 0; a < gnx; a++)
133 for (r = 0; r < DIM; r++)
136 for (c = 0; c < DIM; c++)
138 xf += R[r][c]*xi[j][a][c];
140 msd += w_rls[a]*sqr(xi[i][a][r] - xf);
144 srmsd[i] += sqrt(msd);
145 srmsd[j] += sqrt(msd);
152 min_srmsd = GMX_REAL_MAX;
156 for (i = 0; i < nfr; i++)
158 srmsd[i] /= (nfr - 1);
159 if (srmsd[i] < min_srmsd)
161 min_srmsd = srmsd[i];
165 srmsd_tot += srmsd[i];
169 printf("Average RMSD between all structures: %.3f\n", srmsd_tot/nfr);
170 printf("Structure with lowest RMSD to all others: time %g, av. RMSD %.3f\n",
173 for (a = 0; a < gnx; a++)
175 copy_rvec(xi[min_fr][a], x_ref[index[a]]);
181 int gmx_rotmat(int argc, char *argv[])
183 const char *desc[] = {
184 "[TT]g_rotmat[tt] plots the rotation matrix required for least squares fitting",
185 "a conformation onto the reference conformation provided with",
186 "[TT]-s[tt]. Translation is removed before fitting.",
187 "The output are the three vectors that give the new directions",
188 "of the x, y and z directions of the reference conformation,",
189 "for example: (zx,zy,zz) is the orientation of the reference",
190 "z-axis in the trajectory frame.",
192 "This tool is useful for, for instance,",
193 "determining the orientation of a molecule",
194 "at an interface, possibly on a trajectory produced with",
195 "[TT]trjconv -fit rotxy+transxy[tt] to remove the rotation",
196 "in the [IT]x-y[it] plane.",
198 "Option [TT]-ref[tt] determines a reference structure for fitting,",
199 "instead of using the structure from [TT]-s[tt]. The structure with",
200 "the lowest sum of RMSD's to all other structures is used.",
201 "Since the computational cost of this procedure grows with",
202 "the square of the number of frames, the [TT]-skip[tt] option",
203 "can be useful. A full fit or only a fit in the [IT]x-y[it] plane can",
206 "Option [TT]-fitxy[tt] fits in the [IT]x-y[it] plane before determining",
207 "the rotation matrix."
209 const char *reffit[] =
210 { NULL, "none", "xyz", "xy", NULL };
212 static gmx_bool bFitXY = FALSE, bMW = TRUE;
214 { "-ref", FALSE, etENUM, {reffit},
215 "Determine the optimal reference structure" },
216 { "-skip", FALSE, etINT, {&skip},
217 "Use every nr-th frame for [TT]-ref[tt]" },
218 { "-fitxy", FALSE, etBOOL, {&bFitXY},
219 "Fit the x/y rotation before determining the rotation" },
220 { "-mw", FALSE, etBOOL, {&bMW},
221 "Use mass weighted fitting" }
231 char *grpname, title[256];
233 gmx_rmpbc_t gpbc = NULL;
237 const char *leg[] = { "xx", "xy", "xz", "yx", "yy", "yz", "zx", "zy", "zz" };
238 #define NLEG asize(leg)
240 { efTRX, "-f", NULL, ffREAD },
241 { efTPS, NULL, NULL, ffREAD },
242 { efNDX, NULL, NULL, ffOPTRD },
243 { efXVG, NULL, "rotmat", ffWRITE }
245 #define NFILE asize(fnm)
247 parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW | PCA_BE_NICE,
248 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv);
250 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &x_ref, NULL, box, bMW);
252 gpbc = gmx_rmpbc_init(&top.idef, ePBC, top.atoms.nr, box);
254 gmx_rmpbc(gpbc, top.atoms.nr, box, x_ref);
256 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
258 if (reffit[0][0] != 'n')
260 get_refx(oenv, ftp2fn(efTRX, NFILE, fnm), reffit[0][2] == 'z' ? 3 : 2, skip,
261 gnx, index, bMW, &top, ePBC, x_ref);
264 natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
267 for (i = 0; i < gnx; i++)
269 if (index[i] >= natoms)
271 gmx_fatal(FARGS, "Atom index (%d) is larger than the number of atoms in the trajecory (%d)", index[i]+1, natoms);
273 w_rls[index[i]] = (bMW ? top.atoms.atom[index[i]].m : 1.0);
276 if (reffit[0][0] == 'n')
278 reset_x(gnx, index, natoms, NULL, x_ref, w_rls);
281 out = xvgropen(ftp2fn(efXVG, NFILE, fnm),
282 "Fit matrix", "Time (ps)", "", oenv);
283 xvgr_legend(out, NLEG, leg, oenv);
287 gmx_rmpbc(gpbc, natoms, box, x);
289 reset_x(gnx, index, natoms, NULL, x, w_rls);
293 do_fit_ndim(2, natoms, w_rls, x_ref, x);
296 calc_fit_R(DIM, natoms, w_rls, x_ref, x, R);
299 "%7g %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f\n",
301 R[XX][XX], R[XX][YY], R[XX][ZZ],
302 R[YY][XX], R[YY][YY], R[YY][ZZ],
303 R[ZZ][XX], R[ZZ][YY], R[ZZ][ZZ]);
305 while (read_next_x(oenv, status, &t, natoms, x, box));
307 gmx_rmpbc_done(gpbc);
313 do_view(oenv, ftp2fn(efXVG, NFILE, fnm), "-nxy");