Fixed make_edi.c
[gromacs/rigid-bodies.git] / src / tools / gmx_rotmat.c
bloba0f6d9df20523620a1bd3e71a27388785f90d691
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 bool bMW,t_topology *top,int ePBC,rvec *x_ref)
66 int natoms,nfr_all,nfr,i,j,a,r,c,min_fr;
67 int 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;
75 nfr_all = 0;
76 nfr = 0;
77 snew(ti,100);
78 snew(xi,100);
79 natoms = read_first_x(oenv,&status,trxfn,&ti[nfr],&x,box);
81 snew(w_rls,gnx);
82 tot_mass = 0;
83 for(a=0; a<gnx; a++)
85 if (index[a] >= natoms)
87 gmx_fatal(FARGS,"Atom index (%d) is larger than the number of atoms in the trajecory (%d)",index[a]+1,natoms);
89 w_rls[a] = (bMW ? top->atoms.atom[index[a]].m : 1.0);
90 tot_mass += w_rls[a];
95 if (nfr_all % skip == 0)
97 rm_pbc(&top->idef,ePBC,natoms,box,x,x);
98 snew(xi[nfr],gnx);
99 for(i=0; i<gnx; i++)
101 copy_rvec(x[index[i]],xi[nfr][i]);
103 reset_x(gnx,NULL,gnx,NULL,xi[nfr],w_rls);
104 nfr++;
105 if (nfr % 100 == 0)
107 srenew(ti,nfr+100);
108 srenew(xi,nfr+100);
111 nfr_all++;
113 while(read_next_x(oenv,status,&ti[nfr],natoms,x,box));
114 close_trj(status);
115 sfree(x);
117 snew(srmsd,nfr);
118 for(i=0; i<nfr; i++)
120 printf("\rProcessing frame %d of %d",i,nfr);
121 for(j=i+1; j<nfr; j++)
123 calc_fit_R(nfitdim,gnx,w_rls,xi[i],xi[j],R);
125 msd = 0;
126 for(a=0; a<gnx; a++)
128 for(r=0; r<DIM; r++)
130 xf = 0;
131 for(c=0; c<DIM; c++)
133 xf += R[r][c]*xi[j][a][c];
135 msd += w_rls[a]*sqr(xi[i][a][r] - xf);
138 msd /= tot_mass;
139 srmsd[i] += sqrt(msd);
140 srmsd[j] += sqrt(msd);
142 sfree(xi[i]);
144 printf("\n");
145 sfree(w_rls);
147 min_srmsd = GMX_REAL_MAX;
148 min_fr = -1;
149 min_t = -1;
150 srmsd_tot = 0;
151 for(i=0; i<nfr; i++)
153 srmsd[i] /= (nfr - 1);
154 if (srmsd[i] < min_srmsd)
156 min_srmsd = srmsd[i];
157 min_fr = i;
158 min_t = ti[i];
160 srmsd_tot += srmsd[i];
162 sfree(srmsd);
164 printf("Average RMSD between all structures: %.3f\n",srmsd_tot/nfr);
165 printf("Structure with lowest RMSD to all others: time %g, av. RMSD %.3f\n",
166 min_t,min_srmsd);
168 for(a=0; a<gnx; a++)
170 copy_rvec(xi[min_fr][a],x_ref[index[a]]);
173 sfree(xi);
176 int gmx_rotmat(int argc,char *argv[])
178 const char *desc[] = {
179 "g_rotmat plots the rotation matrix required for least squares fitting",
180 "a conformation onto the reference conformation provided with",
181 "[TT]-s[tt]. Translation is removed before fitting.",
182 "The output are the three vectors that give the new directions",
183 "of the x, y and z directions of the reference conformation,",
184 "for example: (zx,zy,zz) is the orientation of the reference",
185 "z-axis in the trajectory frame.",
186 "[PAR]",
187 "This tool is useful for, for instance,",
188 "determining the orientation of a molecule",
189 "at an interface, possibly on a trajectory produced with",
190 "[TT]trjconv -fit rotxy+transxy[tt] to remove the rotation",
191 "in the xy-plane.",
192 "[PAR]",
193 "Option [TT]-ref[tt] determines a reference structure for fitting,",
194 "instead of using the structure from [TT]-s[tt]. The structure with",
195 "the lowest sum of RMSD's to all other structures is used.",
196 "Since the computational cost of this procedure grows with",
197 "the square of the number of frames, the [TT]-skip[tt] option",
198 "can be useful. A full fit or only a fit in the x/y plane can",
199 "be performed.",
200 "[PAR]",
201 "Option [TT]-fitxy[tt] fits in the x/y plane before determining",
202 "the rotation matrix."
204 const char *reffit[] =
205 { NULL, "none", "xyz", "xy", NULL };
206 static int skip=1;
207 static bool bFitXY=FALSE,bMW=TRUE;
208 t_pargs pa[] = {
209 { "-ref", FALSE, etENUM, {reffit},
210 "Determine the optimal reference structure" },
211 { "-skip", FALSE, etINT, {&skip},
212 "Use every nr-th frame for -ref" },
213 { "-fitxy", FALSE, etBOOL, {&bFitXY},
214 "Fit the x/y rotation before determining the rotation" },
215 { "-mw", FALSE, etBOOL, {&bMW},
216 "Use mass weighted fitting" }
218 FILE *out;
219 int status;
220 t_topology top;
221 int ePBC;
222 rvec *x_ref,*x;
223 matrix box,R;
224 real t;
225 int natoms,i;
226 char *grpname,title[256];
227 int gnx;
228 atom_id *index;
229 output_env_t oenv;
230 real *w_rls;
231 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 }
239 #define NFILE asize(fnm)
241 CopyRight(stderr,argv[0]);
243 parse_common_args(&argc,argv,PCA_CAN_TIME | PCA_CAN_VIEW | PCA_BE_NICE,
244 NFILE,fnm,asize(pa),pa,asize(desc),desc,0,NULL,&oenv);
246 read_tps_conf(ftp2fn(efTPS,NFILE,fnm),title,&top,&ePBC,&x_ref,NULL,box,bMW);
248 rm_pbc(&top.idef,ePBC,top.atoms.nr,box,x_ref,x_ref);
249 get_index(&top.atoms,ftp2fn_null(efNDX,NFILE,fnm),1,&gnx,&index,&grpname);
251 if (reffit[0][0] != 'n')
253 get_refx(oenv,ftp2fn(efTRX,NFILE,fnm),reffit[0][2]=='z' ? 3 : 2,skip,
254 gnx,index,bMW,&top,ePBC,x_ref);
257 natoms = read_first_x(oenv,&status,ftp2fn(efTRX,NFILE,fnm),&t,&x,box);
259 snew(w_rls,natoms);
260 for(i=0; i<gnx; i++)
262 if (index[i] >= natoms)
264 gmx_fatal(FARGS,"Atom index (%d) is larger than the number of atoms in the trajecory (%d)",index[i]+1,natoms);
266 w_rls[index[i]] = (bMW ? top.atoms.atom[index[i]].m : 1.0);
269 if (reffit[0][0] == 'n')
271 reset_x(gnx,index,natoms,NULL,x_ref,w_rls);
274 out = xvgropen(ftp2fn(efXVG,NFILE,fnm),
275 "Fit matrix","Time (ps)","",oenv);
276 xvgr_legend(out,NLEG,leg,oenv);
280 rm_pbc(&top.idef,ePBC,natoms,box,x,x);
282 reset_x(gnx,index,natoms,NULL,x,w_rls);
284 if (bFitXY)
286 do_fit_ndim(2,natoms,w_rls,x_ref,x);
289 calc_fit_R(DIM,natoms,w_rls,x_ref,x,R);
291 fprintf(out,
292 "%7g %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f %7.4f\n",
294 R[XX][XX],R[XX][YY],R[XX][ZZ],
295 R[YY][XX],R[YY][YY],R[YY][ZZ],
296 R[ZZ][XX],R[ZZ][YY],R[ZZ][ZZ]);
298 while(read_next_x(oenv,status,&t,natoms,x,box));
300 close_trj(status);
302 ffclose(out);
304 do_view(oenv,ftp2fn(efXVG,NFILE,fnm),"-nxy");
306 thanx(stderr);
308 return 0;