2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
55 #include "gmx_fatal.h"
61 #include "eigensolver.h"
65 static real find_pdb_bfac(t_atoms *atoms,t_resinfo *ri,char *atomnm)
70 strcpy(rresnm,*ri->name);
72 for(i=0; (i<atoms->nr); i++) {
73 if ((ri->nr == atoms->resinfo[atoms->atom[i].resind].nr) &&
74 (ri->ic == atoms->resinfo[atoms->atom[i].resind].ic) &&
75 (strcmp(*atoms->resinfo[atoms->atom[i].resind].name,rresnm) == 0) &&
76 (strstr(*atoms->atomname[i],atomnm) != NULL))
80 fprintf(stderr,"\rCan not find %s%d-%s in pdbfile\n",
81 rresnm,ri->nr,atomnm);
85 return atoms->pdbinfo[i].bfac;
88 void correlate_aniso(const char *fn,t_atoms *ref,t_atoms *calc,
89 const output_env_t oenv)
94 fp = xvgropen(fn,"Correlation between X-Ray and Computed Uij","X-Ray",
96 for(i=0; (i<ref->nr); i++) {
97 if (ref->pdbinfo[i].bAnisotropic) {
98 for(j=U11; (j<=U23); j++)
99 fprintf(fp,"%10d %10d\n",ref->pdbinfo[i].uij[j],calc->pdbinfo[i].uij[j]);
105 static void average_residues(double f[],double **U,int uind,
106 int isize,atom_id index[],real w_rls[],
115 for(i=0; i<isize; i++) {
116 av += w_rls[index[i]]*(f != NULL ? f[i] : U[i][uind]);
117 m += w_rls[index[i]];
119 atoms->atom[index[i]].resind!=atoms->atom[index[i+1]].resind) {
122 for(j=start; j<=i; j++)
125 for(j=start; j<=i; j++)
135 void print_dir(FILE *fp,real *Uaver)
137 real eigvec[DIM*DIM];
142 fprintf(fp,"MSF X Y Z\n");
145 fprintf(fp," %c ",'X'+d-XX);
147 fprintf(fp," %9.2e",Uaver[3*m+d]);
148 fprintf(fp,"%s\n",m==DIM ? " (nm^2)" : "");
151 for(m=0; m<DIM*DIM; m++)
155 eigensolver(tmp,DIM,0,DIM,eigval,eigvec);
157 fprintf(fp,"\n Eigenvectors\n\n");
158 fprintf(fp,"Eigv %-8.2e %-8.2e %-8.2e (nm^2)\n\n",
159 eigval[2],eigval[1],eigval[0]);
162 fprintf(fp," %c ",'X'+d-XX);
163 for(m=DIM-1; m>=0; m--)
164 fprintf(fp,"%7.4f ",eigvec[3*m+d]);
169 int gmx_rmsf(int argc,char *argv[])
171 const char *desc[] = {
172 "[TT]g_rmsf[tt] computes the root mean square fluctuation (RMSF, i.e. standard ",
173 "deviation) of atomic positions in the trajectory (supplied with [TT]-f[tt])",
174 "after (optionally) fitting to a reference frame (supplied with [TT]-s[tt]).[PAR]",
175 "With option [TT]-oq[tt] the RMSF values are converted to B-factor",
176 "values, which are written to a [TT].pdb[tt] file with the coordinates, of the",
177 "structure file, or of a [TT].pdb[tt] file when [TT]-q[tt] is specified.",
178 "Option [TT]-ox[tt] writes the B-factors to a file with the average",
180 "With the option [TT]-od[tt] the root mean square deviation with",
181 "respect to the reference structure is calculated.[PAR]",
182 "With the option [TT]-aniso[tt], [TT]g_rmsf[tt] will compute anisotropic",
183 "temperature factors and then it will also output average coordinates",
184 "and a [TT].pdb[tt] file with ANISOU records (corresonding to the [TT]-oq[tt]",
185 "or [TT]-ox[tt] option). Please note that the U values",
186 "are orientation-dependent, so before comparison with experimental data",
187 "you should verify that you fit to the experimental coordinates.[PAR]",
188 "When a [TT].pdb[tt] input file is passed to the program and the [TT]-aniso[tt]",
190 "a correlation plot of the Uij will be created, if any anisotropic",
191 "temperature factors are present in the [TT].pdb[tt] file.[PAR]",
192 "With option [TT]-dir[tt] the average MSF (3x3) matrix is diagonalized.",
193 "This shows the directions in which the atoms fluctuate the most and",
196 static gmx_bool bRes=FALSE,bAniso=FALSE,bdevX=FALSE,bFit=TRUE;
198 { "-res", FALSE, etBOOL, {&bRes},
199 "Calculate averages for each residue" },
200 { "-aniso",FALSE, etBOOL, {&bAniso},
201 "Compute anisotropic termperature factors" },
202 { "-fit", FALSE, etBOOL, {&bFit},
203 "Do a least squares superposition before computing RMSF. Without this you must make sure that the reference structure and the trajectory match." }
206 int step,nre,natoms,i,g,m,teller=0;
207 real t,lambda,*w_rls,*w_rms;
213 t_atoms *pdbatoms,*refatoms;
224 FILE *fp; /* the graphics file */
225 const char *devfn,*dirfn;
233 real bfac,pdb_bfac,*Uaver;
237 double *rmsf,invcount,totmass;
241 gmx_rmpbc_t gpbc=NULL;
245 const char *leg[2] = { "MD", "X-Ray" };
248 { efTRX, "-f", NULL, ffREAD },
249 { efTPS, NULL, NULL, ffREAD },
250 { efNDX, NULL, NULL, ffOPTRD },
251 { efPDB, "-q", NULL, ffOPTRD },
252 { efPDB, "-oq", "bfac", ffOPTWR },
253 { efPDB, "-ox", "xaver", ffOPTWR },
254 { efXVG, "-o", "rmsf", ffWRITE },
255 { efXVG, "-od", "rmsdev", ffOPTWR },
256 { efXVG, "-oc", "correl", ffOPTWR },
257 { efLOG, "-dir", "rmsf", ffOPTWR }
259 #define NFILE asize(fnm)
261 CopyRight(stderr,argv[0]);
262 parse_common_args(&argc,argv,PCA_CAN_TIME | PCA_CAN_VIEW | PCA_BE_NICE ,
263 NFILE,fnm,asize(pargs),pargs,asize(desc),desc,0,NULL,
266 bReadPDB = ftp2bSet(efPDB,NFILE,fnm);
267 devfn = opt2fn_null("-od",NFILE,fnm);
268 dirfn = opt2fn_null("-dir",NFILE,fnm);
270 read_tps_conf(ftp2fn(efTPS,NFILE,fnm),title,&top,&ePBC,&xref,NULL,box,TRUE);
271 snew(w_rls,top.atoms.nr);
273 fprintf(stderr,"Select group(s) for root mean square calculation\n");
274 get_index(&top.atoms,ftp2fn_null(efNDX,NFILE,fnm),1,&isize,&index,&grpnames);
277 for(i=0; i<isize; i++)
278 w_rls[index[i]]=top.atoms.atom[index[i]].m;
280 /* Malloc the rmsf arrays */
283 for(i=0; i<isize; i++)
290 get_stx_coordnum(opt2fn("-q",NFILE,fnm),&npdbatoms);
293 init_t_atoms(pdbatoms,npdbatoms,TRUE);
294 init_t_atoms(refatoms,npdbatoms,TRUE);
295 snew(pdbx,npdbatoms);
296 /* Read coordinates twice */
297 read_stx_conf(opt2fn("-q",NFILE,fnm),title,pdbatoms,pdbx,NULL,NULL,pdbbox);
298 read_stx_conf(opt2fn("-q",NFILE,fnm),title,refatoms,pdbx,NULL,NULL,pdbbox);
300 pdbatoms = &top.atoms;
301 refatoms = &top.atoms;
303 npdbatoms = pdbatoms->nr;
304 snew(pdbatoms->pdbinfo,npdbatoms);
305 copy_mat(box,pdbbox);
309 sub_xcm(xref,isize,index,top.atoms.atom,xcm,FALSE);
312 natom = read_first_x(oenv,&status,ftp2fn(efTRX,NFILE,fnm),&t,&x,box);
315 gpbc = gmx_rmpbc_init(&top.idef,ePBC,natom,box);
318 /* Now read the trj again to compute fluctuations */
322 /* Remove periodic boundary */
323 gmx_rmpbc(gpbc,natom,box,x);
325 /* Set center of mass to zero */
326 sub_xcm(x,isize,index,top.atoms.atom,xcm,FALSE);
328 /* Fit to reference structure */
329 do_fit(natom,w_rls,xref,x);
332 /* Calculate Anisotropic U Tensor */
333 for(i=0; i<isize; i++) {
335 for(d=0; d<DIM; d++) {
336 xav[i*DIM + d] += x[aid][d];
338 U[i][d*DIM + m] += x[aid][d]*x[aid][m];
343 /* Calculate RMS Deviation */
344 for(i=0;(i<isize);i++) {
346 for(d=0;(d<DIM);d++) {
347 rmsd_x[i][d] += sqr(x[aid][d]-xref[aid][d]);
353 } while(read_next_x(oenv,status,&t,natom,x,box));
357 gmx_rmpbc_done(gpbc);
360 invcount = 1.0/count;
363 for(i=0; i<isize; i++) {
365 xav[i*DIM + d] *= invcount;
367 for(m=0; m<DIM; m++) {
368 U[i][d*DIM + m] = U[i][d*DIM + m]*invcount
369 - xav[i*DIM + d]*xav[i*DIM + m];
370 Uaver[3*d+m] += top.atoms.atom[index[i]].m*U[i][d*DIM + m];
372 totmass += top.atoms.atom[index[i]].m;
374 for(d=0; d<DIM*DIM; d++)
378 for(d=0; d<DIM*DIM; d++) {
379 average_residues(NULL,U,d,isize,index,w_rls,&top.atoms);
384 for(i=0; i<isize; i++) {
386 pdbatoms->pdbinfo[aid].bAnisotropic = TRUE;
387 pdbatoms->pdbinfo[aid].uij[U11] = 1e6*U[i][XX*DIM + XX];
388 pdbatoms->pdbinfo[aid].uij[U22] = 1e6*U[i][YY*DIM + YY];
389 pdbatoms->pdbinfo[aid].uij[U33] = 1e6*U[i][ZZ*DIM + ZZ];
390 pdbatoms->pdbinfo[aid].uij[U12] = 1e6*U[i][XX*DIM + YY];
391 pdbatoms->pdbinfo[aid].uij[U13] = 1e6*U[i][XX*DIM + ZZ];
392 pdbatoms->pdbinfo[aid].uij[U23] = 1e6*U[i][YY*DIM + ZZ];
400 for(i=0; i<isize; i++)
401 rmsf[i] = U[i][XX*DIM + XX] + U[i][YY*DIM + YY] + U[i][ZZ*DIM + ZZ];
404 fprintf(stdout,"\n");
405 print_dir(stdout,Uaver);
406 fp = ffopen(dirfn,"w");
411 for(i=0; i<isize; i++)
415 /* Write RMSF output */
417 bfac = 8.0*M_PI*M_PI/3.0*100;
418 fp = xvgropen(ftp2fn(efXVG,NFILE,fnm),"B-Factors",
419 label,"(A\\b\\S\\So\\N\\S2\\N)",oenv);
420 xvgr_legend(fp,2,leg,oenv);
421 for(i=0;(i<isize);i++) {
422 if (!bRes || i+1==isize ||
423 top.atoms.atom[index[i]].resind!=top.atoms.atom[index[i+1]].resind) {
424 resind = top.atoms.atom[index[i]].resind;
425 pdb_bfac = find_pdb_bfac(pdbatoms,&top.atoms.resinfo[resind],
426 *(top.atoms.atomname[index[i]]));
428 fprintf(fp,"%5d %10.5f %10.5f\n",
429 bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : i+1,rmsf[i]*bfac,
435 fp = xvgropen(ftp2fn(efXVG,NFILE,fnm),"RMS fluctuation",label,"(nm)",oenv);
436 for(i=0; i<isize; i++)
437 if (!bRes || i+1==isize ||
438 top.atoms.atom[index[i]].resind!=top.atoms.atom[index[i+1]].resind)
439 fprintf(fp,"%5d %8.4f\n",
440 bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : i+1,sqrt(rmsf[i]));
444 for(i=0; i<isize; i++)
445 pdbatoms->pdbinfo[index[i]].bfac = 800*M_PI*M_PI/3.0*rmsf[i];
448 for(i=0; i<isize; i++)
449 rmsf[i] = (rmsd_x[i][XX]+rmsd_x[i][YY]+rmsd_x[i][ZZ])/count;
451 average_residues(rmsf,NULL,0,isize,index,w_rls,&top.atoms);
452 /* Write RMSD output */
453 fp = xvgropen(devfn,"RMS Deviation",label,"(nm)",oenv);
454 for(i=0; i<isize; i++)
455 if (!bRes || i+1==isize ||
456 top.atoms.atom[index[i]].resind!=top.atoms.atom[index[i+1]].resind)
457 fprintf(fp,"%5d %8.4f\n",
458 bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : i+1,sqrt(rmsf[i]));
462 if (opt2bSet("-oq",NFILE,fnm)) {
463 /* Write a .pdb file with B-factors and optionally anisou records */
464 for(i=0; i<isize; i++)
465 rvec_inc(xref[index[i]],xcm);
466 write_sto_conf_indexed(opt2fn("-oq",NFILE,fnm),title,pdbatoms,pdbx,
467 NULL,ePBC,pdbbox,isize,index);
469 if (opt2bSet("-ox",NFILE,fnm)) {
470 /* Misuse xref as a temporary array */
471 for(i=0; i<isize; i++)
473 xref[index[i]][d] = xcm[d] + xav[i*DIM + d];
474 /* Write a .pdb file with B-factors and optionally anisou records */
475 write_sto_conf_indexed(opt2fn("-ox",NFILE,fnm),title,pdbatoms,xref,NULL,
476 ePBC,pdbbox,isize,index);
479 correlate_aniso(opt2fn("-oc",NFILE,fnm),refatoms,pdbatoms,oenv);
480 do_view(oenv,opt2fn("-oc",NFILE,fnm),"-nxy");
482 do_view(oenv,opt2fn("-o",NFILE,fnm),"-nxy");
484 do_view(oenv,opt2fn("-od",NFILE,fnm),"-nxy");