2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
42 #include "gromacs/commandline/pargs.h"
44 #include "gromacs/utility/smalloc.h"
46 #include "gromacs/math/vec.h"
48 #include "gromacs/utility/futil.h"
50 #include "gromacs/fileio/xvgr.h"
52 #include "gromacs/pbcutil/rmpbc.h"
53 #include "gromacs/fileio/tpxio.h"
54 #include "gromacs/fileio/trxio.h"
57 #include "gromacs/math/do_fit.h"
59 static void get_refx(output_env_t oenv, const char *trxfn, int nfitdim, int skip,
61 gmx_bool bMW, t_topology *top, int ePBC, rvec *x_ref)
63 int natoms, nfr_all, nfr, i, j, a, r, c, min_fr;
66 double tot_mass, msd, *srmsd, min_srmsd, srmsd_tot;
71 gmx_rmpbc_t gpbc = NULL;
78 natoms = read_first_x(oenv, &status, trxfn, &ti[nfr], &x, box);
82 for (a = 0; a < gnx; a++)
84 if (index[a] >= natoms)
86 gmx_fatal(FARGS, "Atom index (%d) is larger than the number of atoms in the trajecory (%d)", index[a]+1, natoms);
88 w_rls[a] = (bMW ? top->atoms.atom[index[a]].m : 1.0);
91 gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms);
95 if (nfr_all % skip == 0)
97 gmx_rmpbc(gpbc, natoms, box, x);
99 for (i = 0; i < gnx; i++)
101 copy_rvec(x[index[i]], xi[nfr][i]);
103 reset_x(gnx, NULL, gnx, NULL, xi[nfr], w_rls);
113 while (read_next_x(oenv, status, &ti[nfr], x, box));
117 gmx_rmpbc_done(gpbc);
120 for (i = 0; i < nfr; i++)
122 printf("\rProcessing frame %d of %d", i, nfr);
123 for (j = i+1; j < nfr; j++)
125 calc_fit_R(nfitdim, gnx, w_rls, xi[i], xi[j], R);
128 for (a = 0; a < gnx; a++)
130 for (r = 0; r < DIM; r++)
133 for (c = 0; c < DIM; c++)
135 xf += R[r][c]*xi[j][a][c];
137 msd += w_rls[a]*sqr(xi[i][a][r] - xf);
141 srmsd[i] += sqrt(msd);
142 srmsd[j] += sqrt(msd);
149 min_srmsd = GMX_REAL_MAX;
153 for (i = 0; i < nfr; i++)
155 srmsd[i] /= (nfr - 1);
156 if (srmsd[i] < min_srmsd)
158 min_srmsd = srmsd[i];
162 srmsd_tot += srmsd[i];
166 printf("Average RMSD between all structures: %.3f\n", srmsd_tot/nfr);
167 printf("Structure with lowest RMSD to all others: time %g, av. RMSD %.3f\n",
170 for (a = 0; a < gnx; a++)
172 copy_rvec(xi[min_fr][a], x_ref[index[a]]);
178 int gmx_rotmat(int argc, char *argv[])
180 const char *desc[] = {
181 "[THISMODULE] plots the rotation matrix required for least squares fitting",
182 "a conformation onto the reference conformation provided with",
183 "[TT]-s[tt]. Translation is removed before fitting.",
184 "The output are the three vectors that give the new directions",
185 "of the x, y and z directions of the reference conformation,",
186 "for example: (zx,zy,zz) is the orientation of the reference",
187 "z-axis in the trajectory frame.",
189 "This tool is useful for, for instance,",
190 "determining the orientation of a molecule",
191 "at an interface, possibly on a trajectory produced with",
192 "[TT]gmx trjconv -fit rotxy+transxy[tt] to remove the rotation",
193 "in the [IT]x-y[it] plane.",
195 "Option [TT]-ref[tt] determines a reference structure for fitting,",
196 "instead of using the structure from [TT]-s[tt]. The structure with",
197 "the lowest sum of RMSD's to all other structures is used.",
198 "Since the computational cost of this procedure grows with",
199 "the square of the number of frames, the [TT]-skip[tt] option",
200 "can be useful. A full fit or only a fit in the [IT]x-y[it] plane can",
203 "Option [TT]-fitxy[tt] fits in the [IT]x-y[it] plane before determining",
204 "the rotation matrix."
206 const char *reffit[] =
207 { NULL, "none", "xyz", "xy", NULL };
209 static gmx_bool bFitXY = FALSE, bMW = TRUE;
211 { "-ref", FALSE, etENUM, {reffit},
212 "Determine the optimal reference structure" },
213 { "-skip", FALSE, etINT, {&skip},
214 "Use every nr-th frame for [TT]-ref[tt]" },
215 { "-fitxy", FALSE, etBOOL, {&bFitXY},
216 "Fit the x/y rotation before determining the rotation" },
217 { "-mw", FALSE, etBOOL, {&bMW},
218 "Use mass weighted fitting" }
228 char *grpname, title[256];
230 gmx_rmpbc_t gpbc = NULL;
234 const char *leg[] = { "xx", "xy", "xz", "yx", "yy", "yz", "zx", "zy", "zz" };
235 #define NLEG asize(leg)
237 { efTRX, "-f", NULL, ffREAD },
238 { efTPS, NULL, NULL, ffREAD },
239 { efNDX, NULL, NULL, ffOPTRD },
240 { efXVG, NULL, "rotmat", ffWRITE }
242 #define NFILE asize(fnm)
244 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW | PCA_BE_NICE,
245 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);
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, x, box));
307 gmx_rmpbc_done(gpbc);
313 do_view(oenv, ftp2fn(efXVG, NFILE, fnm), "-nxy");