Merge branch release-4-6
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_rotmat.c
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  *
4  *                This source code is part of
5  *
6  *                 G   R   O   M   A   C   S
7  *
8  *          GROningen MAchine for Chemical Simulations
9  *
10  *                        VERSION 3.2.0
11  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13  * Copyright (c) 2001-2004, The GROMACS development team,
14  * check out http://www.gromacs.org for more information.
15
16  * This program is free software; you can redistribute it and/or
17  * modify it under the terms of the GNU General Public License
18  * as published by the Free Software Foundation; either version 2
19  * of the License, or (at your option) any later version.
20  *
21  * If you want to redistribute modifications, please consider that
22  * scientific software is very special. Version control is crucial -
23  * bugs must be traceable. We will be happy to consider code for
24  * inclusion in the official distribution, but derived work must not
25  * be called official GROMACS. Details are found in the README & COPYING
26  * files - if they are missing, get the official version at www.gromacs.org.
27  *
28  * To help us fund GROMACS development, we humbly ask that you cite
29  * the papers on the package - you can find them in the top README file.
30  *
31  * For more info, check our website at http://www.gromacs.org
32  *
33  * And Hey:
34  * Green Red Orange Magenta Azure Cyan Skyblue
35  */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include <math.h>
41 #include <string.h>
42
43 #include "statutil.h"
44 #include "sysstuff.h"
45 #include "typedefs.h"
46 #include "smalloc.h"
47 #include "macros.h"
48 #include "vec.h"
49 #include "pbc.h"
50 #include "futil.h"
51 #include "statutil.h"
52 #include "index.h"
53 #include "mshift.h"
54 #include "xvgr.h"
55 #include "rmpbc.h"
56 #include "tpxio.h"
57 #include "do_fit.h"
58 #include "gmx_ana.h"
59
60
61 static void get_refx(output_env_t oenv, const char *trxfn, int nfitdim, int skip,
62                      int gnx, int *index,
63                      gmx_bool bMW, t_topology *top, int ePBC, rvec *x_ref)
64 {
65     int          natoms, nfr_all, nfr, i, j, a, r, c, min_fr;
66     t_trxstatus *status;
67     real        *ti, min_t;
68     double       tot_mass, msd, *srmsd, min_srmsd, srmsd_tot;
69     rvec        *x, **xi;
70     real         xf;
71     matrix       box, R;
72     real        *w_rls;
73     gmx_rmpbc_t  gpbc = NULL;
74
75
76     nfr_all = 0;
77     nfr     = 0;
78     snew(ti, 100);
79     snew(xi, 100);
80     natoms = read_first_x(oenv, &status, trxfn, &ti[nfr], &x, box);
81
82     snew(w_rls, gnx);
83     tot_mass = 0;
84     for (a = 0; a < gnx; a++)
85     {
86         if (index[a] >= natoms)
87         {
88             gmx_fatal(FARGS, "Atom index (%d) is larger than the number of atoms in the trajecory (%d)", index[a]+1, natoms);
89         }
90         w_rls[a]  = (bMW ? top->atoms.atom[index[a]].m : 1.0);
91         tot_mass += w_rls[a];
92     }
93     gpbc = gmx_rmpbc_init(&top->idef, ePBC, natoms);
94
95     do
96     {
97         if (nfr_all % skip == 0)
98         {
99             gmx_rmpbc(gpbc, natoms, box, x);
100             snew(xi[nfr], gnx);
101             for (i = 0; i < gnx; i++)
102             {
103                 copy_rvec(x[index[i]], xi[nfr][i]);
104             }
105             reset_x(gnx, NULL, gnx, NULL, xi[nfr], w_rls);
106             nfr++;
107             if (nfr % 100 == 0)
108             {
109                 srenew(ti, nfr+100);
110                 srenew(xi, nfr+100);
111             }
112         }
113         nfr_all++;
114     }
115     while (read_next_x(oenv, status, &ti[nfr], x, box));
116     close_trj(status);
117     sfree(x);
118
119     gmx_rmpbc_done(gpbc);
120
121     snew(srmsd, nfr);
122     for (i = 0; i < nfr; i++)
123     {
124         printf("\rProcessing frame %d of %d", i, nfr);
125         for (j = i+1; j < nfr; j++)
126         {
127             calc_fit_R(nfitdim, gnx, w_rls, xi[i], xi[j], R);
128
129             msd = 0;
130             for (a = 0; a < gnx; a++)
131             {
132                 for (r = 0; r < DIM; r++)
133                 {
134                     xf = 0;
135                     for (c = 0; c < DIM; c++)
136                     {
137                         xf += R[r][c]*xi[j][a][c];
138                     }
139                     msd += w_rls[a]*sqr(xi[i][a][r] - xf);
140                 }
141             }
142             msd      /= tot_mass;
143             srmsd[i] += sqrt(msd);
144             srmsd[j] += sqrt(msd);
145         }
146         sfree(xi[i]);
147     }
148     printf("\n");
149     sfree(w_rls);
150
151     min_srmsd = GMX_REAL_MAX;
152     min_fr    = -1;
153     min_t     = -1;
154     srmsd_tot = 0;
155     for (i = 0; i < nfr; i++)
156     {
157         srmsd[i] /= (nfr - 1);
158         if (srmsd[i] < min_srmsd)
159         {
160             min_srmsd = srmsd[i];
161             min_fr    = i;
162             min_t     = ti[i];
163         }
164         srmsd_tot += srmsd[i];
165     }
166     sfree(srmsd);
167
168     printf("Average RMSD between all structures: %.3f\n", srmsd_tot/nfr);
169     printf("Structure with lowest RMSD to all others: time %g, av. RMSD %.3f\n",
170            min_t, min_srmsd);
171
172     for (a = 0; a < gnx; a++)
173     {
174         copy_rvec(xi[min_fr][a], x_ref[index[a]]);
175     }
176
177     sfree(xi);
178 }
179
180 int gmx_rotmat(int argc, char *argv[])
181 {
182     const char     *desc[] = {
183         "[TT]g_rotmat[tt] plots the rotation matrix required for least squares fitting",
184         "a conformation onto the reference conformation provided with",
185         "[TT]-s[tt]. Translation is removed before fitting.",
186         "The output are the three vectors that give the new directions",
187         "of the x, y and z directions of the reference conformation,",
188         "for example: (zx,zy,zz) is the orientation of the reference",
189         "z-axis in the trajectory frame.",
190         "[PAR]",
191         "This tool is useful for, for instance,",
192         "determining the orientation of a molecule",
193         "at an interface, possibly on a trajectory produced with",
194         "[TT]trjconv -fit rotxy+transxy[tt] to remove the rotation",
195         "in the [IT]x-y[it] plane.",
196         "[PAR]",
197         "Option [TT]-ref[tt] determines a reference structure for fitting,",
198         "instead of using the structure from [TT]-s[tt]. The structure with",
199         "the lowest sum of RMSD's to all other structures is used.",
200         "Since the computational cost of this procedure grows with",
201         "the square of the number of frames, the [TT]-skip[tt] option",
202         "can be useful. A full fit or only a fit in the [IT]x-y[it] plane can",
203         "be performed.",
204         "[PAR]",
205         "Option [TT]-fitxy[tt] fits in the [IT]x-y[it] plane before determining",
206         "the rotation matrix."
207     };
208     const char     *reffit[] =
209     { NULL, "none", "xyz", "xy", NULL };
210     static int      skip   = 1;
211     static gmx_bool bFitXY = FALSE, bMW = TRUE;
212     t_pargs         pa[]   = {
213         { "-ref", FALSE, etENUM, {reffit},
214           "Determine the optimal reference structure" },
215         { "-skip", FALSE, etINT, {&skip},
216           "Use every nr-th frame for [TT]-ref[tt]" },
217         { "-fitxy", FALSE, etBOOL, {&bFitXY},
218           "Fit the x/y rotation before determining the rotation" },
219         { "-mw", FALSE, etBOOL, {&bMW},
220           "Use mass weighted fitting" }
221     };
222     FILE           *out;
223     t_trxstatus    *status;
224     t_topology      top;
225     int             ePBC;
226     rvec           *x_ref, *x;
227     matrix          box, R;
228     real            t;
229     int             natoms, i;
230     char           *grpname, title[256];
231     int             gnx;
232     gmx_rmpbc_t     gpbc = NULL;
233     atom_id        *index;
234     output_env_t    oenv;
235     real           *w_rls;
236     const char     *leg[]  = { "xx", "xy", "xz", "yx", "yy", "yz", "zx", "zy", "zz" };
237 #define NLEG asize(leg)
238     t_filenm        fnm[] = {
239         { efTRX, "-f",   NULL,       ffREAD },
240         { efTPS, NULL,   NULL,       ffREAD },
241         { efNDX, NULL,   NULL,       ffOPTRD },
242         { efXVG, NULL,   "rotmat",   ffWRITE }
243     };
244 #define NFILE asize(fnm)
245
246     if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW | PCA_BE_NICE,
247                            NFILE, fnm, asize(pa), pa, asize(desc), desc, 0, NULL, &oenv))
248     {
249         return 0;
250     }
251
252     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &x_ref, NULL, box, bMW);
253
254     gpbc = gmx_rmpbc_init(&top.idef, ePBC, top.atoms.nr);
255
256     gmx_rmpbc(gpbc, top.atoms.nr, box, x_ref);
257
258     get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &gnx, &index, &grpname);
259
260     if (reffit[0][0] != 'n')
261     {
262         get_refx(oenv, ftp2fn(efTRX, NFILE, fnm), reffit[0][2] == 'z' ? 3 : 2, skip,
263                  gnx, index, bMW, &top, ePBC, x_ref);
264     }
265
266     natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
267
268     snew(w_rls, natoms);
269     for (i = 0; i < gnx; i++)
270     {
271         if (index[i] >= natoms)
272         {
273             gmx_fatal(FARGS, "Atom index (%d) is larger than the number of atoms in the trajecory (%d)", index[i]+1, natoms);
274         }
275         w_rls[index[i]] = (bMW ? top.atoms.atom[index[i]].m : 1.0);
276     }
277
278     if (reffit[0][0] == 'n')
279     {
280         reset_x(gnx, index, natoms, NULL, x_ref, w_rls);
281     }
282
283     out = xvgropen(ftp2fn(efXVG, NFILE, fnm),
284                    "Fit matrix", "Time (ps)", "", oenv);
285     xvgr_legend(out, NLEG, leg, oenv);
286
287     do
288     {
289         gmx_rmpbc(gpbc, natoms, box, x);
290
291         reset_x(gnx, index, natoms, NULL, x, w_rls);
292
293         if (bFitXY)
294         {
295             do_fit_ndim(2, natoms, w_rls, x_ref, x);
296         }
297
298         calc_fit_R(DIM, natoms, w_rls, x_ref, x, R);
299
300         fprintf(out,
301                 "%7g %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f\n",
302                 t,
303                 R[XX][XX], R[XX][YY], R[XX][ZZ],
304                 R[YY][XX], R[YY][YY], R[YY][ZZ],
305                 R[ZZ][XX], R[ZZ][YY], R[ZZ][ZZ]);
306     }
307     while (read_next_x(oenv, status, &t, x, box));
308
309     gmx_rmpbc_done(gpbc);
310
311     close_trj(status);
312
313     ffclose(out);
314
315     do_view(oenv, ftp2fn(efXVG, NFILE, fnm), "-nxy");
316
317     return 0;
318 }