2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009-2017, The GROMACS development team.
5 * Copyright (c) 2019,2020, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
41 #include "gromacs/commandline/pargs.h"
42 #include "gromacs/commandline/viewit.h"
43 #include "gromacs/fileio/confio.h"
44 #include "gromacs/fileio/trxio.h"
45 #include "gromacs/fileio/xvgr.h"
46 #include "gromacs/gmxana/gmx_ana.h"
47 #include "gromacs/math/do_fit.h"
48 #include "gromacs/math/functions.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/pbcutil/rmpbc.h"
51 #include "gromacs/topology/index.h"
52 #include "gromacs/topology/topology.h"
53 #include "gromacs/utility/arraysize.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/futil.h"
56 #include "gromacs/utility/gmxassert.h"
57 #include "gromacs/utility/smalloc.h"
59 static void get_refx(gmx_output_env_t* oenv,
66 const t_topology* top,
70 int natoms, nfr_all, nfr, i, j, a, r, c, min_fr;
73 double tot_mass, msd, *srmsd, min_srmsd, srmsd_tot;
78 gmx_rmpbc_t gpbc = nullptr;
85 natoms = read_first_x(oenv, &status, trxfn, &ti[nfr], &x, box);
89 for (a = 0; a < gnx; a++)
91 if (index[a] >= natoms)
94 "Atom index (%d) is larger than the number of atoms in the trajecory (%d)",
95 index[a] + 1, natoms);
97 w_rls[a] = (bMW ? top->atoms.atom[index[a]].m : 1.0);
100 gpbc = gmx_rmpbc_init(&top->idef, pbcType, natoms);
104 if (nfr_all % skip == 0)
106 gmx_rmpbc(gpbc, natoms, box, x);
108 for (i = 0; i < gnx; i++)
110 copy_rvec(x[index[i]], xi[nfr][i]);
112 reset_x(gnx, nullptr, gnx, nullptr, xi[nfr], w_rls);
116 srenew(ti, nfr + 100);
117 srenew(xi, nfr + 100);
121 } while (read_next_x(oenv, status, &ti[nfr], x, box));
125 gmx_rmpbc_done(gpbc);
128 for (i = 0; i < nfr; i++)
130 fprintf(stdout, "\rProcessing frame %d of %d", i, nfr);
132 for (j = i + 1; j < nfr; j++)
134 calc_fit_R(nfitdim, gnx, w_rls, xi[i], xi[j], R);
137 for (a = 0; a < gnx; a++)
139 for (r = 0; r < DIM; r++)
142 for (c = 0; c < DIM; c++)
144 xf += R[r][c] * xi[j][a][c];
146 msd += w_rls[a] * gmx::square(xi[i][a][r] - xf);
150 srmsd[i] += std::sqrt(msd);
151 srmsd[j] += std::sqrt(msd);
158 min_srmsd = GMX_REAL_MAX;
162 for (i = 0; i < nfr; i++)
164 srmsd[i] /= (nfr - 1);
165 if (srmsd[i] < min_srmsd)
167 min_srmsd = srmsd[i];
171 srmsd_tot += srmsd[i];
175 printf("Average RMSD between all structures: %.3f\n", srmsd_tot / nfr);
176 printf("Structure with lowest RMSD to all others: time %g, av. RMSD %.3f\n", min_t, min_srmsd);
178 for (a = 0; a < gnx; a++)
180 copy_rvec(xi[min_fr][a], x_ref[index[a]]);
186 int gmx_rotmat(int argc, char* argv[])
188 const char* desc[] = {
189 "[THISMODULE] plots the rotation matrix required for least squares fitting",
190 "a conformation onto the reference conformation provided with",
191 "[TT]-s[tt]. Translation is removed before fitting.",
192 "The output are the three vectors that give the new directions",
193 "of the x, y and z directions of the reference conformation,",
194 "for example: (zx,zy,zz) is the orientation of the reference",
195 "z-axis in the trajectory frame.",
197 "This tool is useful for, for instance,",
198 "determining the orientation of a molecule",
199 "at an interface, possibly on a trajectory produced with",
200 "[TT]gmx trjconv -fit rotxy+transxy[tt] to remove the rotation",
201 "in the [IT]x-y[it] plane.",
203 "Option [TT]-ref[tt] determines a reference structure for fitting,",
204 "instead of using the structure from [TT]-s[tt]. The structure with",
205 "the lowest sum of RMSD's to all other structures is used.",
206 "Since the computational cost of this procedure grows with",
207 "the square of the number of frames, the [TT]-skip[tt] option",
208 "can be useful. A full fit or only a fit in the [IT]x-y[it] plane can",
211 "Option [TT]-fitxy[tt] fits in the [IT]x-y[it] plane before determining",
212 "the rotation matrix."
214 const char* reffit[] = { nullptr, "none", "xyz", "xy", nullptr };
216 static gmx_bool bFitXY = FALSE, bMW = TRUE;
218 { "-ref", FALSE, etENUM, { reffit }, "Determine the optimal reference structure" },
219 { "-skip", FALSE, etINT, { &skip }, "Use every nr-th frame for [TT]-ref[tt]" },
224 "Fit the x/y rotation before determining the rotation" },
225 { "-mw", FALSE, etBOOL, { &bMW }, "Use mass weighted fitting" }
237 gmx_rmpbc_t gpbc = nullptr;
239 gmx_output_env_t* oenv;
241 const char* leg[] = { "xx", "xy", "xz", "yx", "yy", "yz", "zx", "zy", "zz" };
242 #define NLEG asize(leg)
243 t_filenm fnm[] = { { efTRX, "-f", nullptr, ffREAD },
244 { efTPS, nullptr, nullptr, ffREAD },
245 { efNDX, nullptr, nullptr, ffOPTRD },
246 { efXVG, nullptr, "rotmat", ffWRITE } };
247 #define NFILE asize(fnm)
249 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW, NFILE, fnm, asize(pa), pa,
250 asize(desc), desc, 0, nullptr, &oenv))
255 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &pbcType, &x_ref, nullptr, box, bMW);
257 gpbc = gmx_rmpbc_init(&top.idef, pbcType, top.atoms.nr);
259 gmx_rmpbc(gpbc, top.atoms.nr, box, x_ref);
261 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
263 GMX_RELEASE_ASSERT(reffit[0] != nullptr, "Options inconsistency; reffit[0] is NULL");
264 if (reffit[0][0] != 'n')
266 get_refx(oenv, ftp2fn(efTRX, NFILE, fnm), reffit[0][2] == 'z' ? 3 : 2, skip, gnx, index,
267 bMW, &top, pbcType, x_ref);
270 natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
273 for (i = 0; i < gnx; i++)
275 if (index[i] >= natoms)
278 "Atom index (%d) is larger than the number of atoms in the trajecory (%d)",
279 index[i] + 1, natoms);
281 w_rls[index[i]] = (bMW ? top.atoms.atom[index[i]].m : 1.0);
284 if (reffit[0][0] == 'n')
286 reset_x(gnx, index, natoms, nullptr, x_ref, w_rls);
289 out = xvgropen(ftp2fn(efXVG, NFILE, fnm), "Fit matrix", "Time (ps)", "", oenv);
290 xvgr_legend(out, NLEG, leg, oenv);
294 gmx_rmpbc(gpbc, natoms, box, x);
296 reset_x(gnx, index, natoms, nullptr, x, w_rls);
300 do_fit_ndim(2, natoms, w_rls, x_ref, x);
303 calc_fit_R(DIM, natoms, w_rls, x_ref, x, R);
305 fprintf(out, "%7g %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f\n", t, R[XX][XX],
306 R[XX][YY], R[XX][ZZ], R[YY][XX], R[YY][YY], R[YY][ZZ], R[ZZ][XX], R[ZZ][YY], R[ZZ][ZZ]);
307 } while (read_next_x(oenv, status, &t, x, box));
309 gmx_rmpbc_done(gpbc);
315 do_view(oenv, ftp2fn(efXVG, NFILE, fnm), "-nxy");