3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * Green Red Orange Magenta Azure Cyan Skyblue
52 #include "gmx_fatal.h"
60 #include "gromacs/linearalgebra/eigensolver.h"
62 static real find_pdb_bfac(t_atoms *atoms, t_resinfo *ri, char *atomnm)
67 strcpy(rresnm, *ri->name);
69 for (i = 0; (i < atoms->nr); i++)
71 if ((ri->nr == atoms->resinfo[atoms->atom[i].resind].nr) &&
72 (ri->ic == atoms->resinfo[atoms->atom[i].resind].ic) &&
73 (strcmp(*atoms->resinfo[atoms->atom[i].resind].name, rresnm) == 0) &&
74 (strstr(*atoms->atomname[i], atomnm) != NULL))
81 fprintf(stderr, "\rCan not find %s%d-%s in pdbfile\n",
82 rresnm, ri->nr, atomnm);
86 return atoms->pdbinfo[i].bfac;
89 void correlate_aniso(const char *fn, t_atoms *ref, t_atoms *calc,
90 const output_env_t oenv)
95 fp = xvgropen(fn, "Correlation between X-Ray and Computed Uij", "X-Ray",
97 for (i = 0; (i < ref->nr); i++)
99 if (ref->pdbinfo[i].bAnisotropic)
101 for (j = U11; (j <= U23); j++)
103 fprintf(fp, "%10d %10d\n", ref->pdbinfo[i].uij[j], calc->pdbinfo[i].uij[j]);
110 static void average_residues(double f[], double **U, int uind,
111 int isize, atom_id index[], real w_rls[],
120 for (i = 0; i < isize; i++)
122 av += w_rls[index[i]]*(f != NULL ? f[i] : U[i][uind]);
123 m += w_rls[index[i]];
125 atoms->atom[index[i]].resind != atoms->atom[index[i+1]].resind)
130 for (j = start; j <= i; j++)
137 for (j = start; j <= i; j++)
149 void print_dir(FILE *fp, real *Uaver)
151 real eigvec[DIM*DIM];
156 fprintf(fp, "MSF X Y Z\n");
157 for (d = 0; d < DIM; d++)
159 fprintf(fp, " %c ", 'X'+d-XX);
160 for (m = 0; m < DIM; m++)
162 fprintf(fp, " %9.2e", Uaver[3*m+d]);
164 fprintf(fp, "%s\n", m == DIM ? " (nm^2)" : "");
167 for (m = 0; m < DIM*DIM; m++)
173 eigensolver(tmp, DIM, 0, DIM, eigval, eigvec);
175 fprintf(fp, "\n Eigenvectors\n\n");
176 fprintf(fp, "Eigv %-8.2e %-8.2e %-8.2e (nm^2)\n\n",
177 eigval[2], eigval[1], eigval[0]);
178 for (d = 0; d < DIM; d++)
180 fprintf(fp, " %c ", 'X'+d-XX);
181 for (m = DIM-1; m >= 0; m--)
183 fprintf(fp, "%7.4f ", eigvec[3*m+d]);
189 int gmx_rmsf(int argc, char *argv[])
191 const char *desc[] = {
192 "[TT]g_rmsf[tt] computes the root mean square fluctuation (RMSF, i.e. standard ",
193 "deviation) of atomic positions in the trajectory (supplied with [TT]-f[tt])",
194 "after (optionally) fitting to a reference frame (supplied with [TT]-s[tt]).[PAR]",
195 "With option [TT]-oq[tt] the RMSF values are converted to B-factor",
196 "values, which are written to a [TT].pdb[tt] file with the coordinates, of the",
197 "structure file, or of a [TT].pdb[tt] file when [TT]-q[tt] is specified.",
198 "Option [TT]-ox[tt] writes the B-factors to a file with the average",
200 "With the option [TT]-od[tt] the root mean square deviation with",
201 "respect to the reference structure is calculated.[PAR]",
202 "With the option [TT]-aniso[tt], [TT]g_rmsf[tt] will compute anisotropic",
203 "temperature factors and then it will also output average coordinates",
204 "and a [TT].pdb[tt] file with ANISOU records (corresonding to the [TT]-oq[tt]",
205 "or [TT]-ox[tt] option). Please note that the U values",
206 "are orientation-dependent, so before comparison with experimental data",
207 "you should verify that you fit to the experimental coordinates.[PAR]",
208 "When a [TT].pdb[tt] input file is passed to the program and the [TT]-aniso[tt]",
210 "a correlation plot of the Uij will be created, if any anisotropic",
211 "temperature factors are present in the [TT].pdb[tt] file.[PAR]",
212 "With option [TT]-dir[tt] the average MSF (3x3) matrix is diagonalized.",
213 "This shows the directions in which the atoms fluctuate the most and",
216 static gmx_bool bRes = FALSE, bAniso = FALSE, bdevX = FALSE, bFit = TRUE;
218 { "-res", FALSE, etBOOL, {&bRes},
219 "Calculate averages for each residue" },
220 { "-aniso", FALSE, etBOOL, {&bAniso},
221 "Compute anisotropic termperature factors" },
222 { "-fit", FALSE, etBOOL, {&bFit},
223 "Do a least squares superposition before computing RMSF. Without this you must make sure that the reference structure and the trajectory match." }
226 int step, nre, natoms, i, g, m, teller = 0;
227 real t, lambda, *w_rls, *w_rms;
233 t_atoms *pdbatoms, *refatoms;
237 rvec *x, *pdbx, *xref;
244 FILE *fp; /* the graphics file */
245 const char *devfn, *dirfn;
253 real bfac, pdb_bfac, *Uaver;
257 double *rmsf, invcount, totmass;
261 gmx_rmpbc_t gpbc = NULL;
265 const char *leg[2] = { "MD", "X-Ray" };
268 { efTRX, "-f", NULL, ffREAD },
269 { efTPS, NULL, NULL, ffREAD },
270 { efNDX, NULL, NULL, ffOPTRD },
271 { efPDB, "-q", NULL, ffOPTRD },
272 { efPDB, "-oq", "bfac", ffOPTWR },
273 { efPDB, "-ox", "xaver", ffOPTWR },
274 { efXVG, "-o", "rmsf", ffWRITE },
275 { efXVG, "-od", "rmsdev", ffOPTWR },
276 { efXVG, "-oc", "correl", ffOPTWR },
277 { efLOG, "-dir", "rmsf", ffOPTWR }
279 #define NFILE asize(fnm)
281 parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW | PCA_BE_NICE,
282 NFILE, fnm, asize(pargs), pargs, asize(desc), desc, 0, NULL,
285 bReadPDB = ftp2bSet(efPDB, NFILE, fnm);
286 devfn = opt2fn_null("-od", NFILE, fnm);
287 dirfn = opt2fn_null("-dir", NFILE, fnm);
289 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xref, NULL, box, TRUE);
290 snew(w_rls, top.atoms.nr);
292 fprintf(stderr, "Select group(s) for root mean square calculation\n");
293 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpnames);
296 for (i = 0; i < isize; i++)
298 w_rls[index[i]] = top.atoms.atom[index[i]].m;
301 /* Malloc the rmsf arrays */
302 snew(xav, isize*DIM);
304 for (i = 0; i < isize; i++)
316 get_stx_coordnum(opt2fn("-q", NFILE, fnm), &npdbatoms);
319 init_t_atoms(pdbatoms, npdbatoms, TRUE);
320 init_t_atoms(refatoms, npdbatoms, TRUE);
321 snew(pdbx, npdbatoms);
322 /* Read coordinates twice */
323 read_stx_conf(opt2fn("-q", NFILE, fnm), title, pdbatoms, pdbx, NULL, NULL, pdbbox);
324 read_stx_conf(opt2fn("-q", NFILE, fnm), title, refatoms, pdbx, NULL, NULL, pdbbox);
328 pdbatoms = &top.atoms;
329 refatoms = &top.atoms;
331 npdbatoms = pdbatoms->nr;
332 snew(pdbatoms->pdbinfo, npdbatoms);
333 copy_mat(box, pdbbox);
338 sub_xcm(xref, isize, index, top.atoms.atom, xcm, FALSE);
341 natom = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
345 gpbc = gmx_rmpbc_init(&top.idef, ePBC, natom, box);
348 /* Now read the trj again to compute fluctuations */
354 /* Remove periodic boundary */
355 gmx_rmpbc(gpbc, natom, box, x);
357 /* Set center of mass to zero */
358 sub_xcm(x, isize, index, top.atoms.atom, xcm, FALSE);
360 /* Fit to reference structure */
361 do_fit(natom, w_rls, xref, x);
364 /* Calculate Anisotropic U Tensor */
365 for (i = 0; i < isize; i++)
368 for (d = 0; d < DIM; d++)
370 xav[i*DIM + d] += x[aid][d];
371 for (m = 0; m < DIM; m++)
373 U[i][d*DIM + m] += x[aid][d]*x[aid][m];
380 /* Calculate RMS Deviation */
381 for (i = 0; (i < isize); i++)
384 for (d = 0; (d < DIM); d++)
386 rmsd_x[i][d] += sqr(x[aid][d]-xref[aid][d]);
393 while (read_next_x(oenv, status, &t, natom, x, box));
398 gmx_rmpbc_done(gpbc);
402 invcount = 1.0/count;
403 snew(Uaver, DIM*DIM);
405 for (i = 0; i < isize; i++)
407 for (d = 0; d < DIM; d++)
409 xav[i*DIM + d] *= invcount;
411 for (d = 0; d < DIM; d++)
413 for (m = 0; m < DIM; m++)
415 U[i][d*DIM + m] = U[i][d*DIM + m]*invcount
416 - xav[i*DIM + d]*xav[i*DIM + m];
417 Uaver[3*d+m] += top.atoms.atom[index[i]].m*U[i][d*DIM + m];
420 totmass += top.atoms.atom[index[i]].m;
422 for (d = 0; d < DIM*DIM; d++)
429 for (d = 0; d < DIM*DIM; d++)
431 average_residues(NULL, U, d, isize, index, w_rls, &top.atoms);
437 for (i = 0; i < isize; i++)
440 pdbatoms->pdbinfo[aid].bAnisotropic = TRUE;
441 pdbatoms->pdbinfo[aid].uij[U11] = 1e6*U[i][XX*DIM + XX];
442 pdbatoms->pdbinfo[aid].uij[U22] = 1e6*U[i][YY*DIM + YY];
443 pdbatoms->pdbinfo[aid].uij[U33] = 1e6*U[i][ZZ*DIM + ZZ];
444 pdbatoms->pdbinfo[aid].uij[U12] = 1e6*U[i][XX*DIM + YY];
445 pdbatoms->pdbinfo[aid].uij[U13] = 1e6*U[i][XX*DIM + ZZ];
446 pdbatoms->pdbinfo[aid].uij[U23] = 1e6*U[i][YY*DIM + ZZ];
458 for (i = 0; i < isize; i++)
460 rmsf[i] = U[i][XX*DIM + XX] + U[i][YY*DIM + YY] + U[i][ZZ*DIM + ZZ];
465 fprintf(stdout, "\n");
466 print_dir(stdout, Uaver);
467 fp = ffopen(dirfn, "w");
468 print_dir(fp, Uaver);
472 for (i = 0; i < isize; i++)
478 /* Write RMSF output */
481 bfac = 8.0*M_PI*M_PI/3.0*100;
482 fp = xvgropen(ftp2fn(efXVG, NFILE, fnm), "B-Factors",
483 label, "(A\\b\\S\\So\\N\\S2\\N)", oenv);
484 xvgr_legend(fp, 2, leg, oenv);
485 for (i = 0; (i < isize); i++)
487 if (!bRes || i+1 == isize ||
488 top.atoms.atom[index[i]].resind != top.atoms.atom[index[i+1]].resind)
490 resind = top.atoms.atom[index[i]].resind;
491 pdb_bfac = find_pdb_bfac(pdbatoms, &top.atoms.resinfo[resind],
492 *(top.atoms.atomname[index[i]]));
494 fprintf(fp, "%5d %10.5f %10.5f\n",
495 bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : index[i]+1, rmsf[i]*bfac,
503 fp = xvgropen(ftp2fn(efXVG, NFILE, fnm), "RMS fluctuation", label, "(nm)", oenv);
504 for (i = 0; i < isize; i++)
506 if (!bRes || i+1 == isize ||
507 top.atoms.atom[index[i]].resind != top.atoms.atom[index[i+1]].resind)
509 fprintf(fp, "%5d %8.4f\n",
510 bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : index[i]+1, sqrt(rmsf[i]));
516 for (i = 0; i < isize; i++)
518 pdbatoms->pdbinfo[index[i]].bfac = 800*M_PI*M_PI/3.0*rmsf[i];
523 for (i = 0; i < isize; i++)
525 rmsf[i] = (rmsd_x[i][XX]+rmsd_x[i][YY]+rmsd_x[i][ZZ])/count;
529 average_residues(rmsf, NULL, 0, isize, index, w_rls, &top.atoms);
531 /* Write RMSD output */
532 fp = xvgropen(devfn, "RMS Deviation", label, "(nm)", oenv);
533 for (i = 0; i < isize; i++)
535 if (!bRes || i+1 == isize ||
536 top.atoms.atom[index[i]].resind != top.atoms.atom[index[i+1]].resind)
538 fprintf(fp, "%5d %8.4f\n",
539 bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : index[i]+1, sqrt(rmsf[i]));
545 if (opt2bSet("-oq", NFILE, fnm))
547 /* Write a .pdb file with B-factors and optionally anisou records */
548 for (i = 0; i < isize; i++)
550 rvec_inc(xref[index[i]], xcm);
552 write_sto_conf_indexed(opt2fn("-oq", NFILE, fnm), title, pdbatoms, pdbx,
553 NULL, ePBC, pdbbox, isize, index);
555 if (opt2bSet("-ox", NFILE, fnm))
557 /* Misuse xref as a temporary array */
558 for (i = 0; i < isize; i++)
560 for (d = 0; d < DIM; d++)
562 xref[index[i]][d] = xcm[d] + xav[i*DIM + d];
565 /* Write a .pdb file with B-factors and optionally anisou records */
566 write_sto_conf_indexed(opt2fn("-ox", NFILE, fnm), title, pdbatoms, xref, NULL,
567 ePBC, pdbbox, isize, index);
571 correlate_aniso(opt2fn("-oc", NFILE, fnm), refatoms, pdbatoms, oenv);
572 do_view(oenv, opt2fn("-oc", NFILE, fnm), "-nxy");
574 do_view(oenv, opt2fn("-o", NFILE, fnm), "-nxy");
577 do_view(oenv, opt2fn("-od", NFILE, fnm), "-nxy");