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.
40 #include "gromacs/commandline/pargs.h"
42 #include "gromacs/utility/smalloc.h"
44 #include "gromacs/math/vec.h"
45 #include "gromacs/utility/futil.h"
46 #include "gromacs/topology/index.h"
47 #include "gromacs/fileio/xvgr.h"
49 #include "gromacs/pbcutil/rmpbc.h"
50 #include "gromacs/fileio/tpxio.h"
51 #include "gromacs/fileio/trxio.h"
54 #include "gromacs/math/do_fit.h"
56 static void get_refx(output_env_t oenv, const char *trxfn, int nfitdim, int skip,
58 gmx_bool bMW, t_topology *top, int ePBC, rvec *x_ref)
60 int natoms, nfr_all, nfr, i, j, a, r, c, min_fr;
63 double tot_mass, msd, *srmsd, min_srmsd, srmsd_tot;
68 gmx_rmpbc_t gpbc = NULL;
75 natoms = read_first_x(oenv, &status, trxfn, &ti[nfr], &x, box);
79 for (a = 0; a < gnx; a++)
81 if (index[a] >= natoms)
83 gmx_fatal(FARGS, "Atom index (%d) is larger than the number of atoms in the trajecory (%d)", index[a]+1, natoms);
85 w_rls[a] = (bMW ? top->atoms.atom[index[a]].m : 1.0);
88 gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms);
92 if (nfr_all % skip == 0)
94 gmx_rmpbc(gpbc, natoms, box, x);
96 for (i = 0; i < gnx; i++)
98 copy_rvec(x[index[i]], xi[nfr][i]);
100 reset_x(gnx, NULL, gnx, NULL, xi[nfr], w_rls);
110 while (read_next_x(oenv, status, &ti[nfr], x, box));
114 gmx_rmpbc_done(gpbc);
117 for (i = 0; i < nfr; i++)
119 printf("\rProcessing frame %d of %d", i, nfr);
120 for (j = i+1; j < nfr; j++)
122 calc_fit_R(nfitdim, gnx, w_rls, xi[i], xi[j], R);
125 for (a = 0; a < gnx; a++)
127 for (r = 0; r < DIM; r++)
130 for (c = 0; c < DIM; c++)
132 xf += R[r][c]*xi[j][a][c];
134 msd += w_rls[a]*sqr(xi[i][a][r] - xf);
138 srmsd[i] += sqrt(msd);
139 srmsd[j] += sqrt(msd);
146 min_srmsd = GMX_REAL_MAX;
150 for (i = 0; i < nfr; i++)
152 srmsd[i] /= (nfr - 1);
153 if (srmsd[i] < min_srmsd)
155 min_srmsd = srmsd[i];
159 srmsd_tot += srmsd[i];
163 printf("Average RMSD between all structures: %.3f\n", srmsd_tot/nfr);
164 printf("Structure with lowest RMSD to all others: time %g, av. RMSD %.3f\n",
167 for (a = 0; a < gnx; a++)
169 copy_rvec(xi[min_fr][a], x_ref[index[a]]);
175 int gmx_rotmat(int argc, char *argv[])
177 const char *desc[] = {
178 "[THISMODULE] plots the rotation matrix required for least squares fitting",
179 "a conformation onto the reference conformation provided with",
180 "[TT]-s[tt]. Translation is removed before fitting.",
181 "The output are the three vectors that give the new directions",
182 "of the x, y and z directions of the reference conformation,",
183 "for example: (zx,zy,zz) is the orientation of the reference",
184 "z-axis in the trajectory frame.",
186 "This tool is useful for, for instance,",
187 "determining the orientation of a molecule",
188 "at an interface, possibly on a trajectory produced with",
189 "[TT]gmx trjconv -fit rotxy+transxy[tt] to remove the rotation",
190 "in the [IT]x-y[it] plane.",
192 "Option [TT]-ref[tt] determines a reference structure for fitting,",
193 "instead of using the structure from [TT]-s[tt]. The structure with",
194 "the lowest sum of RMSD's to all other structures is used.",
195 "Since the computational cost of this procedure grows with",
196 "the square of the number of frames, the [TT]-skip[tt] option",
197 "can be useful. A full fit or only a fit in the [IT]x-y[it] plane can",
200 "Option [TT]-fitxy[tt] fits in the [IT]x-y[it] plane before determining",
201 "the rotation matrix."
203 const char *reffit[] =
204 { NULL, "none", "xyz", "xy", NULL };
206 static gmx_bool bFitXY = FALSE, bMW = TRUE;
208 { "-ref", FALSE, etENUM, {reffit},
209 "Determine the optimal reference structure" },
210 { "-skip", FALSE, etINT, {&skip},
211 "Use every nr-th frame for [TT]-ref[tt]" },
212 { "-fitxy", FALSE, etBOOL, {&bFitXY},
213 "Fit the x/y rotation before determining the rotation" },
214 { "-mw", FALSE, etBOOL, {&bMW},
215 "Use mass weighted fitting" }
225 char *grpname, title[256];
227 gmx_rmpbc_t gpbc = NULL;
231 const char *leg[] = { "xx", "xy", "xz", "yx", "yy", "yz", "zx", "zy", "zz" };
232 #define NLEG asize(leg)
234 { efTRX, "-f", NULL, ffREAD },
235 { efTPS, NULL, NULL, ffREAD },
236 { efNDX, NULL, NULL, ffOPTRD },
237 { efXVG, NULL, "rotmat", ffWRITE }
239 #define NFILE asize(fnm)
241 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW | PCA_BE_NICE,
242 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
247 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &x_ref, NULL, box, bMW);
249 gpbc = gmx_rmpbc_init(&top.idef, ePBC, top.atoms.nr);
251 gmx_rmpbc(gpbc, top.atoms.nr, box, x_ref);
253 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
255 if (reffit[0][0] != 'n')
257 get_refx(oenv, ftp2fn(efTRX, NFILE, fnm), reffit[0][2] == 'z' ? 3 : 2, skip,
258 gnx, index, bMW, &top, ePBC, x_ref);
261 natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
264 for (i = 0; i < gnx; i++)
266 if (index[i] >= natoms)
268 gmx_fatal(FARGS, "Atom index (%d) is larger than the number of atoms in the trajecory (%d)", index[i]+1, natoms);
270 w_rls[index[i]] = (bMW ? top.atoms.atom[index[i]].m : 1.0);
273 if (reffit[0][0] == 'n')
275 reset_x(gnx, index, natoms, NULL, x_ref, w_rls);
278 out = xvgropen(ftp2fn(efXVG, NFILE, fnm),
279 "Fit matrix", "Time (ps)", "", oenv);
280 xvgr_legend(out, NLEG, leg, oenv);
284 gmx_rmpbc(gpbc, natoms, box, x);
286 reset_x(gnx, index, natoms, NULL, x, w_rls);
290 do_fit_ndim(2, natoms, w_rls, x_ref, x);
293 calc_fit_R(DIM, natoms, w_rls, x_ref, x, R);
296 "%7g %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f\n",
298 R[XX][XX], R[XX][YY], R[XX][ZZ],
299 R[YY][XX], R[YY][YY], R[YY][ZZ],
300 R[ZZ][XX], R[ZZ][YY], R[ZZ][ZZ]);
302 while (read_next_x(oenv, status, &t, x, box));
304 gmx_rmpbc_done(gpbc);
310 do_view(oenv, ftp2fn(efXVG, NFILE, fnm), "-nxy");