Cleaned up memory usage in gmx traj and trjconv
[gromacs.git] / src / gromacs / gmxana / gmx_rotmat.cpp
blobea5c2f4222625b627893aa0a8d99d8d4c05b7dc6
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014,2015,2016,2017, 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.
35 #include "gmxpre.h"
37 #include <cmath>
38 #include <cstring>
40 #include "gromacs/commandline/pargs.h"
41 #include "gromacs/commandline/viewit.h"
42 #include "gromacs/fileio/confio.h"
43 #include "gromacs/fileio/trxio.h"
44 #include "gromacs/fileio/xvgr.h"
45 #include "gromacs/gmxana/gmx_ana.h"
46 #include "gromacs/math/do_fit.h"
47 #include "gromacs/math/functions.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/pbcutil/rmpbc.h"
50 #include "gromacs/topology/index.h"
51 #include "gromacs/topology/topology.h"
52 #include "gromacs/utility/arraysize.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/futil.h"
55 #include "gromacs/utility/gmxassert.h"
56 #include "gromacs/utility/smalloc.h"
58 static void get_refx(gmx_output_env_t *oenv, const char *trxfn, int nfitdim, int skip,
59 int gnx, int *index,
60 gmx_bool bMW, const t_topology *top, int ePBC, rvec *x_ref)
62 int natoms, nfr_all, nfr, i, j, a, r, c, min_fr;
63 t_trxstatus *status;
64 real *ti, min_t;
65 double tot_mass, msd, *srmsd, min_srmsd, srmsd_tot;
66 rvec *x, **xi;
67 real xf;
68 matrix box, R;
69 real *w_rls;
70 gmx_rmpbc_t gpbc = nullptr;
73 nfr_all = 0;
74 nfr = 0;
75 snew(ti, 100);
76 snew(xi, 100);
77 natoms = read_first_x(oenv, &status, trxfn, &ti[nfr], &x, box);
79 snew(w_rls, gnx);
80 tot_mass = 0;
81 for (a = 0; a < gnx; a++)
83 if (index[a] >= natoms)
85 gmx_fatal(FARGS, "Atom index (%d) is larger than the number of atoms in the trajecory (%d)", index[a]+1, natoms);
87 w_rls[a] = (bMW ? top->atoms.atom[index[a]].m : 1.0);
88 tot_mass += w_rls[a];
90 gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms);
94 if (nfr_all % skip == 0)
96 gmx_rmpbc(gpbc, natoms, box, x);
97 snew(xi[nfr], gnx);
98 for (i = 0; i < gnx; i++)
100 copy_rvec(x[index[i]], xi[nfr][i]);
102 reset_x(gnx, nullptr, gnx, nullptr, xi[nfr], w_rls);
103 nfr++;
104 if (nfr % 100 == 0)
106 srenew(ti, nfr+100);
107 srenew(xi, nfr+100);
110 nfr_all++;
112 while (read_next_x(oenv, status, &ti[nfr], x, box));
113 close_trx(status);
114 sfree(x);
116 gmx_rmpbc_done(gpbc);
118 snew(srmsd, nfr);
119 for (i = 0; i < nfr; i++)
121 fprintf(stdout, "\rProcessing frame %d of %d", i, nfr);
122 fflush(stdout);
123 for (j = i+1; j < nfr; j++)
125 calc_fit_R(nfitdim, gnx, w_rls, xi[i], xi[j], R);
127 msd = 0;
128 for (a = 0; a < gnx; a++)
130 for (r = 0; r < DIM; r++)
132 xf = 0;
133 for (c = 0; c < DIM; c++)
135 xf += R[r][c]*xi[j][a][c];
137 msd += w_rls[a]*gmx::square(xi[i][a][r] - xf);
140 msd /= tot_mass;
141 srmsd[i] += std::sqrt(msd);
142 srmsd[j] += std::sqrt(msd);
144 sfree(xi[i]);
146 printf("\n");
147 sfree(w_rls);
149 min_srmsd = GMX_REAL_MAX;
150 min_fr = -1;
151 min_t = -1;
152 srmsd_tot = 0;
153 for (i = 0; i < nfr; i++)
155 srmsd[i] /= (nfr - 1);
156 if (srmsd[i] < min_srmsd)
158 min_srmsd = srmsd[i];
159 min_fr = i;
160 min_t = ti[i];
162 srmsd_tot += srmsd[i];
164 sfree(srmsd);
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",
168 min_t, min_srmsd);
170 for (a = 0; a < gnx; a++)
172 copy_rvec(xi[min_fr][a], x_ref[index[a]]);
175 sfree(xi);
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.",
188 "[PAR]",
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.",
194 "[PAR]",
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",
201 "be performed.",
202 "[PAR]",
203 "Option [TT]-fitxy[tt] fits in the [IT]x-y[it] plane before determining",
204 "the rotation matrix."
206 const char *reffit[] =
207 { nullptr, "none", "xyz", "xy", nullptr };
208 static int skip = 1;
209 static gmx_bool bFitXY = FALSE, bMW = TRUE;
210 t_pargs pa[] = {
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" }
220 FILE *out;
221 t_trxstatus *status;
222 t_topology top;
223 int ePBC;
224 rvec *x_ref, *x;
225 matrix box, R;
226 real t;
227 int natoms, i;
228 char *grpname;
229 int gnx;
230 gmx_rmpbc_t gpbc = nullptr;
231 int *index;
232 gmx_output_env_t *oenv;
233 real *w_rls;
234 const char *leg[] = { "xx", "xy", "xz", "yx", "yy", "yz", "zx", "zy", "zz" };
235 #define NLEG asize(leg)
236 t_filenm fnm[] = {
237 { efTRX, "-f", nullptr, ffREAD },
238 { efTPS, nullptr, nullptr, ffREAD },
239 { efNDX, nullptr, nullptr, ffOPTRD },
240 { efXVG, nullptr, "rotmat", ffWRITE }
242 #define NFILE asize(fnm)
244 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW,
245 NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, nullptr, &oenv))
247 return 0;
250 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &ePBC, &x_ref, nullptr, 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 GMX_RELEASE_ASSERT(reffit[0] != nullptr, "Options inconsistency; reffit[0] is NULL");
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);
267 snew(w_rls, natoms);
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, nullptr, 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, nullptr, x, w_rls);
292 if (bFitXY)
294 do_fit_ndim(2, natoms, w_rls, x_ref, x);
297 calc_fit_R(DIM, natoms, w_rls, x_ref, x, R);
299 fprintf(out,
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);
310 close_trx(status);
312 xvgrclose(out);
314 do_view(oenv, ftp2fn(efXVG, NFILE, fnm), "-nxy");
316 return 0;