Redefine the default boolean type to gmx_bool.
[alexxy/gromacs.git] / src / tools / 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 "copyrite.h"
51 #include "futil.h"
52 #include "statutil.h"
53 #include "index.h"
54 #include "mshift.h"
55 #include "xvgr.h"
56 #include "rmpbc.h"
57 #include "tpxio.h"
58 #include "do_fit.h"
59 #include "gmx_ana.h"
60
61
62 static void get_refx(output_env_t oenv,const char *trxfn,int nfitdim,int skip,
63                      int gnx,int *index,
64                      gmx_bool bMW,t_topology *top,int ePBC,rvec *x_ref)
65 {
66     int    natoms,nfr_all,nfr,i,j,a,r,c,min_fr;
67     t_trxstatus *status;
68     real   *ti,min_t;
69     double tot_mass,msd,*srmsd,min_srmsd,srmsd_tot;
70     rvec   *x,**xi;
71     real   xf;
72     matrix box,R;
73     real   *w_rls;
74     gmx_rmpbc_t  gpbc=NULL;
75
76
77     nfr_all = 0;
78     nfr     = 0;
79     snew(ti,100);
80     snew(xi,100);
81     natoms = read_first_x(oenv,&status,trxfn,&ti[nfr],&x,box); 
82
83     snew(w_rls,gnx);
84     tot_mass = 0;
85     for(a=0; a<gnx; a++)
86     {
87         if (index[a] >= natoms)
88         {
89             gmx_fatal(FARGS,"Atom index (%d) is larger than the number of atoms in the trajecory (%d)",index[a]+1,natoms);
90         }
91         w_rls[a] = (bMW ? top->atoms.atom[index[a]].m : 1.0);
92         tot_mass += w_rls[a];
93     }
94     gpbc = gmx_rmpbc_init(&top->idef,ePBC,natoms,box);
95
96     do
97     {
98         if (nfr_all % skip == 0)
99         {
100             gmx_rmpbc(gpbc,natoms,box,x);
101             snew(xi[nfr],gnx);
102             for(i=0; i<gnx; i++)
103             {
104                 copy_rvec(x[index[i]],xi[nfr][i]);
105             }
106             reset_x(gnx,NULL,gnx,NULL,xi[nfr],w_rls);
107             nfr++;
108             if (nfr % 100 == 0)
109             {
110                 srenew(ti,nfr+100);
111                 srenew(xi,nfr+100);
112             }
113         }
114         nfr_all++;
115     }
116     while(read_next_x(oenv,status,&ti[nfr],natoms,x,box));
117     close_trj(status);
118     sfree(x);
119
120     gmx_rmpbc_done(gpbc);
121
122     snew(srmsd,nfr);
123     for(i=0; i<nfr; i++)
124     {
125         printf("\rProcessing frame %d of %d",i,nfr);
126         for(j=i+1; j<nfr; j++)
127         {
128             calc_fit_R(nfitdim,gnx,w_rls,xi[i],xi[j],R);
129
130             msd = 0;
131             for(a=0; a<gnx; a++)
132             {
133                 for(r=0; r<DIM; r++)
134                 {
135                     xf = 0;
136                     for(c=0; c<DIM; c++)
137                     {
138                         xf += R[r][c]*xi[j][a][c];
139                     }
140                     msd += w_rls[a]*sqr(xi[i][a][r] - xf);
141                 }
142             }
143             msd /= tot_mass;
144             srmsd[i] += sqrt(msd);
145             srmsd[j] += sqrt(msd);
146         }
147         sfree(xi[i]);
148     }
149     printf("\n");
150     sfree(w_rls);
151
152     min_srmsd = GMX_REAL_MAX;
153     min_fr    = -1;
154     min_t     = -1;
155     srmsd_tot = 0;
156     for(i=0; i<nfr; i++)
157     {
158         srmsd[i] /= (nfr - 1);
159         if (srmsd[i] < min_srmsd)
160         {
161             min_srmsd = srmsd[i];
162             min_fr    = i;
163             min_t     = ti[i];
164         }
165         srmsd_tot += srmsd[i];
166     }
167     sfree(srmsd);
168
169     printf("Average RMSD between all structures: %.3f\n",srmsd_tot/nfr);
170     printf("Structure with lowest RMSD to all others: time %g, av. RMSD %.3f\n",
171            min_t,min_srmsd);
172
173     for(a=0; a<gnx; a++)
174     {
175         copy_rvec(xi[min_fr][a],x_ref[index[a]]);
176     }
177
178     sfree(xi);
179 }
180
181 int gmx_rotmat(int argc,char *argv[])
182 {
183     const char *desc[] = {
184         "g_rotmat plots the rotation matrix required for least squares fitting",
185         "a conformation onto the reference conformation provided with",
186         "[TT]-s[tt]. Translation is removed before fitting.",
187         "The output are the three vectors that give the new directions",
188         "of the x, y and z directions of the reference conformation,",
189         "for example: (zx,zy,zz) is the orientation of the reference",
190         "z-axis in the trajectory frame.",
191         "[PAR]",
192         "This tool is useful for, for instance,",
193         "determining the orientation of a molecule",
194         "at an interface, possibly on a trajectory produced with",
195         "[TT]trjconv -fit rotxy+transxy[tt] to remove the rotation",
196         "in the xy-plane.",
197         "[PAR]",
198         "Option [TT]-ref[tt] determines a reference structure for fitting,",
199         "instead of using the structure from [TT]-s[tt]. The structure with",
200         "the lowest sum of RMSD's to all other structures is used.",
201         "Since the computational cost of this procedure grows with",
202         "the square of the number of frames, the [TT]-skip[tt] option",
203         "can be useful. A full fit or only a fit in the x/y plane can",
204         "be performed.",
205         "[PAR]",
206         "Option [TT]-fitxy[tt] fits in the x/y plane before determining",
207         "the rotation matrix."
208     };
209     const char *reffit[] = 
210         { NULL, "none", "xyz", "xy", NULL }; 
211     static int  skip=1;
212     static gmx_bool bFitXY=FALSE,bMW=TRUE;
213     t_pargs pa[] = {
214         { "-ref", FALSE, etENUM, {reffit},
215           "Determine the optimal reference structure" },
216         { "-skip", FALSE, etINT, {&skip},
217           "Use every nr-th frame for -ref" },
218         { "-fitxy", FALSE, etBOOL, {&bFitXY},
219           "Fit the x/y rotation before determining the rotation" },
220         { "-mw", FALSE, etBOOL, {&bMW},
221           "Use mass weighted fitting" }
222     };
223     FILE       *out;
224     t_trxstatus *status;
225     t_topology top;
226     int        ePBC;
227     rvec       *x_ref,*x;
228     matrix     box,R;
229     real       t;
230     int        natoms,i;
231     char       *grpname,title[256];
232     int        gnx;
233     gmx_rmpbc_t  gpbc=NULL;
234     atom_id    *index;
235     output_env_t oenv;
236     real       *w_rls;
237     const char *leg[]  = { "xx", "xy", "xz", "yx", "yy", "yz", "zx", "zy", "zz" }; 
238 #define NLEG asize(leg) 
239     t_filenm fnm[] = {
240         { efTRX, "-f",   NULL,       ffREAD }, 
241         { efTPS, NULL,   NULL,       ffREAD },
242         { efNDX, NULL,   NULL,       ffOPTRD },
243         { efXVG, NULL,   "rotmat",   ffWRITE }
244     }; 
245 #define NFILE asize(fnm) 
246     
247     CopyRight(stderr,argv[0]);
248     
249     parse_common_args(&argc,argv,PCA_CAN_TIME | PCA_CAN_VIEW | PCA_BE_NICE,
250                       NFILE,fnm,asize(pa),pa,asize(desc),desc,0,NULL,&oenv); 
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,box);
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,natoms,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     thanx(stderr);
318     
319     return 0;
320 }