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"
45 #include "gromacs/utility/smalloc.h"
49 #include "gromacs/fileio/futil.h"
54 #include "gromacs/fileio/tpxio.h"
55 #include "gromacs/fileio/trxio.h"
58 #include "gromacs/math/do_fit.h"
60 static void get_refx(output_env_t oenv, const char *trxfn, int nfitdim, int skip,
62 gmx_bool bMW, t_topology *top, int ePBC, rvec *x_ref)
64 int natoms, nfr_all, nfr, i, j, a, r, c, min_fr;
67 double tot_mass, msd, *srmsd, min_srmsd, srmsd_tot;
72 gmx_rmpbc_t gpbc = NULL;
79 natoms = read_first_x(oenv, &status, trxfn, &ti[nfr], &x, box);
83 for (a = 0; a < gnx; a++)
85 if (index[a] >= natoms)
87 gmx_fatal(FARGS, "Atom index (%d) is larger than the number of atoms in the trajecory (%d)", index[a]+1, natoms);
89 w_rls[a] = (bMW ? top->atoms.atom[index[a]].m : 1.0);
92 gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms);
96 if (nfr_all % skip == 0)
98 gmx_rmpbc(gpbc, natoms, box, x);
100 for (i = 0; i < gnx; i++)
102 copy_rvec(x[index[i]], xi[nfr][i]);
104 reset_x(gnx, NULL, gnx, NULL, xi[nfr], w_rls);
114 while (read_next_x(oenv, status, &ti[nfr], x, box));
118 gmx_rmpbc_done(gpbc);
121 for (i = 0; i < nfr; i++)
123 printf("\rProcessing frame %d of %d", i, nfr);
124 for (j = i+1; j < nfr; j++)
126 calc_fit_R(nfitdim, gnx, w_rls, xi[i], xi[j], R);
129 for (a = 0; a < gnx; a++)
131 for (r = 0; r < DIM; r++)
134 for (c = 0; c < DIM; c++)
136 xf += R[r][c]*xi[j][a][c];
138 msd += w_rls[a]*sqr(xi[i][a][r] - xf);
142 srmsd[i] += sqrt(msd);
143 srmsd[j] += sqrt(msd);
150 min_srmsd = GMX_REAL_MAX;
154 for (i = 0; i < nfr; i++)
156 srmsd[i] /= (nfr - 1);
157 if (srmsd[i] < min_srmsd)
159 min_srmsd = srmsd[i];
163 srmsd_tot += srmsd[i];
167 printf("Average RMSD between all structures: %.3f\n", srmsd_tot/nfr);
168 printf("Structure with lowest RMSD to all others: time %g, av. RMSD %.3f\n",
171 for (a = 0; a < gnx; a++)
173 copy_rvec(xi[min_fr][a], x_ref[index[a]]);
179 int gmx_rotmat(int argc, char *argv[])
181 const char *desc[] = {
182 "[THISMODULE] plots the rotation matrix required for least squares fitting",
183 "a conformation onto the reference conformation provided with",
184 "[TT]-s[tt]. Translation is removed before fitting.",
185 "The output are the three vectors that give the new directions",
186 "of the x, y and z directions of the reference conformation,",
187 "for example: (zx,zy,zz) is the orientation of the reference",
188 "z-axis in the trajectory frame.",
190 "This tool is useful for, for instance,",
191 "determining the orientation of a molecule",
192 "at an interface, possibly on a trajectory produced with",
193 "[TT]gmx trjconv -fit rotxy+transxy[tt] to remove the rotation",
194 "in the [IT]x-y[it] plane.",
196 "Option [TT]-ref[tt] determines a reference structure for fitting,",
197 "instead of using the structure from [TT]-s[tt]. The structure with",
198 "the lowest sum of RMSD's to all other structures is used.",
199 "Since the computational cost of this procedure grows with",
200 "the square of the number of frames, the [TT]-skip[tt] option",
201 "can be useful. A full fit or only a fit in the [IT]x-y[it] plane can",
204 "Option [TT]-fitxy[tt] fits in the [IT]x-y[it] plane before determining",
205 "the rotation matrix."
207 const char *reffit[] =
208 { NULL, "none", "xyz", "xy", NULL };
210 static gmx_bool bFitXY = FALSE, bMW = TRUE;
212 { "-ref", FALSE, etENUM, {reffit},
213 "Determine the optimal reference structure" },
214 { "-skip", FALSE, etINT, {&skip},
215 "Use every nr-th frame for [TT]-ref[tt]" },
216 { "-fitxy", FALSE, etBOOL, {&bFitXY},
217 "Fit the x/y rotation before determining the rotation" },
218 { "-mw", FALSE, etBOOL, {&bMW},
219 "Use mass weighted fitting" }
229 char *grpname, title[256];
231 gmx_rmpbc_t gpbc = NULL;
235 const char *leg[] = { "xx", "xy", "xz", "yx", "yy", "yz", "zx", "zy", "zz" };
236 #define NLEG asize(leg)
238 { efTRX, "-f", NULL, ffREAD },
239 { efTPS, NULL, NULL, ffREAD },
240 { efNDX, NULL, NULL, ffOPTRD },
241 { efXVG, NULL, "rotmat", ffWRITE }
243 #define NFILE asize(fnm)
245 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW | PCA_BE_NICE,
246 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
251 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &x_ref, NULL, box, bMW);
253 gpbc = gmx_rmpbc_init(&top.idef, ePBC, top.atoms.nr);
255 gmx_rmpbc(gpbc, top.atoms.nr, box, x_ref);
257 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
259 if (reffit[0][0] != 'n')
261 get_refx(oenv, ftp2fn(efTRX, NFILE, fnm), reffit[0][2] == 'z' ? 3 : 2, skip,
262 gnx, index, bMW, &top, ePBC, x_ref);
265 natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
268 for (i = 0; i < gnx; i++)
270 if (index[i] >= natoms)
272 gmx_fatal(FARGS, "Atom index (%d) is larger than the number of atoms in the trajecory (%d)", index[i]+1, natoms);
274 w_rls[index[i]] = (bMW ? top.atoms.atom[index[i]].m : 1.0);
277 if (reffit[0][0] == 'n')
279 reset_x(gnx, index, natoms, NULL, x_ref, w_rls);
282 out = xvgropen(ftp2fn(efXVG, NFILE, fnm),
283 "Fit matrix", "Time (ps)", "", oenv);
284 xvgr_legend(out, NLEG, leg, oenv);
288 gmx_rmpbc(gpbc, natoms, box, x);
290 reset_x(gnx, index, natoms, NULL, x, w_rls);
294 do_fit_ndim(2, natoms, w_rls, x_ref, x);
297 calc_fit_R(DIM, natoms, w_rls, x_ref, x, R);
300 "%7g %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f\n",
302 R[XX][XX], R[XX][YY], R[XX][ZZ],
303 R[YY][XX], R[YY][YY], R[YY][ZZ],
304 R[ZZ][XX], R[ZZ][YY], R[ZZ][ZZ]);
306 while (read_next_x(oenv, status, &t, x, box));
308 gmx_rmpbc_done(gpbc);
314 do_view(oenv, ftp2fn(efXVG, NFILE, fnm), "-nxy");