Merge branch 'master' of git://git.gromacs.org/gromacs
[gromacs/adressmacs.git] / src / tools / gmx_rotmat.c
blob68ffcc6d0ac3ee1bb14b97be43381d1f802b8e64
1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
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.
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.
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.
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.
31 * For more info, check our website at http://www.gromacs.org
33 * And Hey:
34 * Green Red Orange Magenta Azure Cyan Skyblue
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
40 #include <math.h>
41 #include <string.h>
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"
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)
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;
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);
83 snew(w_rls,gnx);
84 tot_mass = 0;
85 for(a=0; a<gnx; a++)
87 if (index[a] >= natoms)
89 gmx_fatal(FARGS,"Atom index (%d) is larger than the number of atoms in the trajecory (%d)",index[a]+1,natoms);
91 w_rls[a] = (bMW ? top->atoms.atom[index[a]].m : 1.0);
92 tot_mass += w_rls[a];
94 gpbc = gmx_rmpbc_init(&top->idef,ePBC,natoms,box);
98 if (nfr_all % skip == 0)
100 gmx_rmpbc(gpbc,natoms,box,x);
101 snew(xi[nfr],gnx);
102 for(i=0; i<gnx; i++)
104 copy_rvec(x[index[i]],xi[nfr][i]);
106 reset_x(gnx,NULL,gnx,NULL,xi[nfr],w_rls);
107 nfr++;
108 if (nfr % 100 == 0)
110 srenew(ti,nfr+100);
111 srenew(xi,nfr+100);
114 nfr_all++;
116 while(read_next_x(oenv,status,&ti[nfr],natoms,x,box));
117 close_trj(status);
118 sfree(x);
120 gmx_rmpbc_done(gpbc);
122 snew(srmsd,nfr);
123 for(i=0; i<nfr; i++)
125 printf("\rProcessing frame %d of %d",i,nfr);
126 for(j=i+1; j<nfr; j++)
128 calc_fit_R(nfitdim,gnx,w_rls,xi[i],xi[j],R);
130 msd = 0;
131 for(a=0; a<gnx; a++)
133 for(r=0; r<DIM; r++)
135 xf = 0;
136 for(c=0; c<DIM; c++)
138 xf += R[r][c]*xi[j][a][c];
140 msd += w_rls[a]*sqr(xi[i][a][r] - xf);
143 msd /= tot_mass;
144 srmsd[i] += sqrt(msd);
145 srmsd[j] += sqrt(msd);
147 sfree(xi[i]);
149 printf("\n");
150 sfree(w_rls);
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++)
158 srmsd[i] /= (nfr - 1);
159 if (srmsd[i] < min_srmsd)
161 min_srmsd = srmsd[i];
162 min_fr = i;
163 min_t = ti[i];
165 srmsd_tot += srmsd[i];
167 sfree(srmsd);
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);
173 for(a=0; a<gnx; a++)
175 copy_rvec(xi[min_fr][a],x_ref[index[a]]);
178 sfree(xi);
181 int gmx_rotmat(int argc,char *argv[])
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."
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" }
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 }
245 #define NFILE asize(fnm)
247 CopyRight(stderr,argv[0]);
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);
252 read_tps_conf(ftp2fn(efTPS,NFILE,fnm),title,&top,&ePBC,&x_ref,NULL,box,bMW);
254 gpbc = gmx_rmpbc_init(&top.idef,ePBC,top.atoms.nr,box);
256 gmx_rmpbc(gpbc,top.atoms.nr,box,x_ref);
258 get_index(&top.atoms,ftp2fn_null(efNDX,NFILE,fnm),1,&gnx,&index,&grpname);
260 if (reffit[0][0] != 'n')
262 get_refx(oenv,ftp2fn(efTRX,NFILE,fnm),reffit[0][2]=='z' ? 3 : 2,skip,
263 gnx,index,bMW,&top,ePBC,x_ref);
266 natoms = read_first_x(oenv,&status,ftp2fn(efTRX,NFILE,fnm),&t,&x,box);
268 snew(w_rls,natoms);
269 for(i=0; i<gnx; i++)
271 if (index[i] >= natoms)
273 gmx_fatal(FARGS,"Atom index (%d) is larger than the number of atoms in the trajecory (%d)",index[i]+1,natoms);
275 w_rls[index[i]] = (bMW ? top.atoms.atom[index[i]].m : 1.0);
278 if (reffit[0][0] == 'n')
280 reset_x(gnx,index,natoms,NULL,x_ref,w_rls);
283 out = xvgropen(ftp2fn(efXVG,NFILE,fnm),
284 "Fit matrix","Time (ps)","",oenv);
285 xvgr_legend(out,NLEG,leg,oenv);
289 gmx_rmpbc(gpbc,natoms,box,x);
291 reset_x(gnx,index,natoms,NULL,x,w_rls);
293 if (bFitXY)
295 do_fit_ndim(2,natoms,w_rls,x_ref,x);
298 calc_fit_R(DIM,natoms,w_rls,x_ref,x,R);
300 fprintf(out,
301 "%7g %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f\n",
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]);
307 while(read_next_x(oenv,status,&t,natoms,x,box));
309 gmx_rmpbc_done(gpbc);
311 close_trj(status);
313 ffclose(out);
315 do_view(oenv,ftp2fn(efXVG,NFILE,fnm),"-nxy");
317 thanx(stderr);
319 return 0;