e4606039b4fc5433773f42297e8d0f825f9c4447
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_rotmat.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 #include "gmxpre.h"
36
37 #include <math.h>
38 #include <string.h>
39
40 #include "gromacs/commandline/pargs.h"
41 #include "gromacs/legacyheaders/typedefs.h"
42 #include "gromacs/utility/smalloc.h"
43 #include "gromacs/legacyheaders/macros.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"
48 #include "gromacs/legacyheaders/viewit.h"
49 #include "gromacs/pbcutil/rmpbc.h"
50 #include "gromacs/fileio/tpxio.h"
51 #include "gromacs/fileio/trxio.h"
52 #include "gmx_ana.h"
53
54 #include "gromacs/math/do_fit.h"
55
56 static void get_refx(output_env_t oenv, const char *trxfn, int nfitdim, int skip,
57                      int gnx, int *index,
58                      gmx_bool bMW, t_topology *top, int ePBC, rvec *x_ref)
59 {
60     int          natoms, nfr_all, nfr, i, j, a, r, c, min_fr;
61     t_trxstatus *status;
62     real        *ti, min_t;
63     double       tot_mass, msd, *srmsd, min_srmsd, srmsd_tot;
64     rvec        *x, **xi;
65     real         xf;
66     matrix       box, R;
67     real        *w_rls;
68     gmx_rmpbc_t  gpbc = NULL;
69
70
71     nfr_all = 0;
72     nfr     = 0;
73     snew(ti, 100);
74     snew(xi, 100);
75     natoms = read_first_x(oenv, &status, trxfn, &ti[nfr], &x, box);
76
77     snew(w_rls, gnx);
78     tot_mass = 0;
79     for (a = 0; a < gnx; a++)
80     {
81         if (index[a] >= natoms)
82         {
83             gmx_fatal(FARGS, "Atom index (%d) is larger than the number of atoms in the trajecory (%d)", index[a]+1, natoms);
84         }
85         w_rls[a]  = (bMW ? top->atoms.atom[index[a]].m : 1.0);
86         tot_mass += w_rls[a];
87     }
88     gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms);
89
90     do
91     {
92         if (nfr_all % skip == 0)
93         {
94             gmx_rmpbc(gpbc, natoms, box, x);
95             snew(xi[nfr], gnx);
96             for (i = 0; i < gnx; i++)
97             {
98                 copy_rvec(x[index[i]], xi[nfr][i]);
99             }
100             reset_x(gnx, NULL, gnx, NULL, xi[nfr], w_rls);
101             nfr++;
102             if (nfr % 100 == 0)
103             {
104                 srenew(ti, nfr+100);
105                 srenew(xi, nfr+100);
106             }
107         }
108         nfr_all++;
109     }
110     while (read_next_x(oenv, status, &ti[nfr], x, box));
111     close_trj(status);
112     sfree(x);
113
114     gmx_rmpbc_done(gpbc);
115
116     snew(srmsd, nfr);
117     for (i = 0; i < nfr; i++)
118     {
119         printf("\rProcessing frame %d of %d", i, nfr);
120         for (j = i+1; j < nfr; j++)
121         {
122             calc_fit_R(nfitdim, gnx, w_rls, xi[i], xi[j], R);
123
124             msd = 0;
125             for (a = 0; a < gnx; a++)
126             {
127                 for (r = 0; r < DIM; r++)
128                 {
129                     xf = 0;
130                     for (c = 0; c < DIM; c++)
131                     {
132                         xf += R[r][c]*xi[j][a][c];
133                     }
134                     msd += w_rls[a]*sqr(xi[i][a][r] - xf);
135                 }
136             }
137             msd      /= tot_mass;
138             srmsd[i] += sqrt(msd);
139             srmsd[j] += sqrt(msd);
140         }
141         sfree(xi[i]);
142     }
143     printf("\n");
144     sfree(w_rls);
145
146     min_srmsd = GMX_REAL_MAX;
147     min_fr    = -1;
148     min_t     = -1;
149     srmsd_tot = 0;
150     for (i = 0; i < nfr; i++)
151     {
152         srmsd[i] /= (nfr - 1);
153         if (srmsd[i] < min_srmsd)
154         {
155             min_srmsd = srmsd[i];
156             min_fr    = i;
157             min_t     = ti[i];
158         }
159         srmsd_tot += srmsd[i];
160     }
161     sfree(srmsd);
162
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",
165            min_t, min_srmsd);
166
167     for (a = 0; a < gnx; a++)
168     {
169         copy_rvec(xi[min_fr][a], x_ref[index[a]]);
170     }
171
172     sfree(xi);
173 }
174
175 int gmx_rotmat(int argc, char *argv[])
176 {
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.",
185         "[PAR]",
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.",
191         "[PAR]",
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",
198         "be performed.",
199         "[PAR]",
200         "Option [TT]-fitxy[tt] fits in the [IT]x-y[it] plane before determining",
201         "the rotation matrix."
202     };
203     const char     *reffit[] =
204     { NULL, "none", "xyz", "xy", NULL };
205     static int      skip   = 1;
206     static gmx_bool bFitXY = FALSE, bMW = TRUE;
207     t_pargs         pa[]   = {
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" }
216     };
217     FILE           *out;
218     t_trxstatus    *status;
219     t_topology      top;
220     int             ePBC;
221     rvec           *x_ref, *x;
222     matrix          box, R;
223     real            t;
224     int             natoms, i;
225     char           *grpname, title[256];
226     int             gnx;
227     gmx_rmpbc_t     gpbc = NULL;
228     atom_id        *index;
229     output_env_t    oenv;
230     real           *w_rls;
231     const char     *leg[]  = { "xx", "xy", "xz", "yx", "yy", "yz", "zx", "zy", "zz" };
232 #define NLEG asize(leg)
233     t_filenm        fnm[] = {
234         { efTRX, "-f",   NULL,       ffREAD },
235         { efTPS, NULL,   NULL,       ffREAD },
236         { efNDX, NULL,   NULL,       ffOPTRD },
237         { efXVG, NULL,   "rotmat",   ffWRITE }
238     };
239 #define NFILE asize(fnm)
240
241     if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW,
242                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
243     {
244         return 0;
245     }
246
247     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &x_ref, NULL, box, bMW);
248
249     gpbc = gmx_rmpbc_init(&top.idef, ePBC, top.atoms.nr);
250
251     gmx_rmpbc(gpbc, top.atoms.nr, box, x_ref);
252
253     get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
254
255     if (reffit[0][0] != 'n')
256     {
257         get_refx(oenv, ftp2fn(efTRX, NFILE, fnm), reffit[0][2] == 'z' ? 3 : 2, skip,
258                  gnx, index, bMW, &top, ePBC, x_ref);
259     }
260
261     natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
262
263     snew(w_rls, natoms);
264     for (i = 0; i < gnx; i++)
265     {
266         if (index[i] >= natoms)
267         {
268             gmx_fatal(FARGS, "Atom index (%d) is larger than the number of atoms in the trajecory (%d)", index[i]+1, natoms);
269         }
270         w_rls[index[i]] = (bMW ? top.atoms.atom[index[i]].m : 1.0);
271     }
272
273     if (reffit[0][0] == 'n')
274     {
275         reset_x(gnx, index, natoms, NULL, x_ref, w_rls);
276     }
277
278     out = xvgropen(ftp2fn(efXVG, NFILE, fnm),
279                    "Fit matrix", "Time (ps)", "", oenv);
280     xvgr_legend(out, NLEG, leg, oenv);
281
282     do
283     {
284         gmx_rmpbc(gpbc, natoms, box, x);
285
286         reset_x(gnx, index, natoms, NULL, x, w_rls);
287
288         if (bFitXY)
289         {
290             do_fit_ndim(2, natoms, w_rls, x_ref, x);
291         }
292
293         calc_fit_R(DIM, natoms, w_rls, x_ref, x, R);
294
295         fprintf(out,
296                 "%7g %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f\n",
297                 t,
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]);
301     }
302     while (read_next_x(oenv, status, &t, x, box));
303
304     gmx_rmpbc_done(gpbc);
305
306     close_trj(status);
307
308     gmx_ffclose(out);
309
310     do_view(oenv, ftp2fn(efXVG, NFILE, fnm), "-nxy");
311
312     return 0;
313 }