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 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
41 #include "gromacs/commandline/pargs.h"
42 #include "gromacs/fileio/confio.h"
43 #include "gromacs/fileio/pdbio.h"
44 #include "gromacs/fileio/tpxio.h"
45 #include "gromacs/fileio/trxio.h"
46 #include "gromacs/fileio/xvgr.h"
47 #include "gromacs/gmxana/gmx_ana.h"
48 #include "gromacs/gmxana/princ.h"
49 #include "gromacs/legacyheaders/macros.h"
50 #include "gromacs/legacyheaders/typedefs.h"
51 #include "gromacs/legacyheaders/viewit.h"
52 #include "gromacs/linearalgebra/eigensolver.h"
53 #include "gromacs/math/do_fit.h"
54 #include "gromacs/math/vec.h"
55 #include "gromacs/pbcutil/rmpbc.h"
56 #include "gromacs/topology/index.h"
57 #include "gromacs/utility/futil.h"
58 #include "gromacs/utility/smalloc.h"
60 static real find_pdb_bfac(t_atoms *atoms, t_resinfo *ri, char *atomnm)
65 strcpy(rresnm, *ri->name);
67 for (i = 0; (i < atoms->nr); i++)
69 if ((ri->nr == atoms->resinfo[atoms->atom[i].resind].nr) &&
70 (ri->ic == atoms->resinfo[atoms->atom[i].resind].ic) &&
71 (strcmp(*atoms->resinfo[atoms->atom[i].resind].name, rresnm) == 0) &&
72 (strstr(*atoms->atomname[i], atomnm) != NULL))
79 fprintf(stderr, "\rCan not find %s%d-%s in pdbfile\n",
80 rresnm, ri->nr, atomnm);
84 return atoms->pdbinfo[i].bfac;
87 void correlate_aniso(const char *fn, t_atoms *ref, t_atoms *calc,
88 const output_env_t oenv)
93 fp = xvgropen(fn, "Correlation between X-Ray and Computed Uij", "X-Ray",
95 for (i = 0; (i < ref->nr); i++)
97 if (ref->pdbinfo[i].bAnisotropic)
99 for (j = U11; (j <= U23); j++)
101 fprintf(fp, "%10d %10d\n", ref->pdbinfo[i].uij[j], calc->pdbinfo[i].uij[j]);
108 static void average_residues(double f[], double **U, int uind,
109 int isize, atom_id index[], real w_rls[],
118 for (i = 0; i < isize; i++)
120 av += w_rls[index[i]]*(f != NULL ? f[i] : U[i][uind]);
121 m += w_rls[index[i]];
123 atoms->atom[index[i]].resind != atoms->atom[index[i+1]].resind)
128 for (j = start; j <= i; j++)
135 for (j = start; j <= i; j++)
147 void print_dir(FILE *fp, real *Uaver)
149 real eigvec[DIM*DIM];
154 fprintf(fp, "MSF X Y Z\n");
155 for (d = 0; d < DIM; d++)
157 fprintf(fp, " %c ", 'X'+d-XX);
158 for (m = 0; m < DIM; m++)
160 fprintf(fp, " %9.2e", Uaver[3*m+d]);
162 fprintf(fp, "%s\n", m == DIM ? " (nm^2)" : "");
165 for (m = 0; m < DIM*DIM; m++)
171 eigensolver(tmp, DIM, 0, DIM, eigval, eigvec);
173 fprintf(fp, "\n Eigenvectors\n\n");
174 fprintf(fp, "Eigv %-8.2e %-8.2e %-8.2e (nm^2)\n\n",
175 eigval[2], eigval[1], eigval[0]);
176 for (d = 0; d < DIM; d++)
178 fprintf(fp, " %c ", 'X'+d-XX);
179 for (m = DIM-1; m >= 0; m--)
181 fprintf(fp, "%7.4f ", eigvec[3*m+d]);
187 int gmx_rmsf(int argc, char *argv[])
189 const char *desc[] = {
190 "[THISMODULE] computes the root mean square fluctuation (RMSF, i.e. standard ",
191 "deviation) of atomic positions in the trajectory (supplied with [TT]-f[tt])",
192 "after (optionally) fitting to a reference frame (supplied with [TT]-s[tt]).[PAR]",
193 "With option [TT]-oq[tt] the RMSF values are converted to B-factor",
194 "values, which are written to a [TT].pdb[tt] file with the coordinates, of the",
195 "structure file, or of a [TT].pdb[tt] file when [TT]-q[tt] is specified.",
196 "Option [TT]-ox[tt] writes the B-factors to a file with the average",
198 "With the option [TT]-od[tt] the root mean square deviation with",
199 "respect to the reference structure is calculated.[PAR]",
200 "With the option [TT]-aniso[tt], [THISMODULE] will compute anisotropic",
201 "temperature factors and then it will also output average coordinates",
202 "and a [TT].pdb[tt] file with ANISOU records (corresonding to the [TT]-oq[tt]",
203 "or [TT]-ox[tt] option). Please note that the U values",
204 "are orientation-dependent, so before comparison with experimental data",
205 "you should verify that you fit to the experimental coordinates.[PAR]",
206 "When a [TT].pdb[tt] input file is passed to the program and the [TT]-aniso[tt]",
208 "a correlation plot of the Uij will be created, if any anisotropic",
209 "temperature factors are present in the [TT].pdb[tt] file.[PAR]",
210 "With option [TT]-dir[tt] the average MSF (3x3) matrix is diagonalized.",
211 "This shows the directions in which the atoms fluctuate the most and",
214 static gmx_bool bRes = FALSE, bAniso = FALSE, bdevX = FALSE, bFit = TRUE;
216 { "-res", FALSE, etBOOL, {&bRes},
217 "Calculate averages for each residue" },
218 { "-aniso", FALSE, etBOOL, {&bAniso},
219 "Compute anisotropic termperature factors" },
220 { "-fit", FALSE, etBOOL, {&bFit},
221 "Do a least squares superposition before computing RMSF. Without this you must make sure that the reference structure and the trajectory match." }
224 int step, nre, natoms, i, g, m, teller = 0;
225 real t, lambda, *w_rls, *w_rms;
231 t_atoms *pdbatoms, *refatoms;
235 rvec *x, *pdbx, *xref;
242 FILE *fp; /* the graphics file */
243 const char *devfn, *dirfn;
251 real bfac, pdb_bfac, *Uaver;
255 double *rmsf, invcount, totmass;
259 gmx_rmpbc_t gpbc = NULL;
263 const char *leg[2] = { "MD", "X-Ray" };
266 { efTRX, "-f", NULL, ffREAD },
267 { efTPS, NULL, NULL, ffREAD },
268 { efNDX, NULL, NULL, ffOPTRD },
269 { efPDB, "-q", NULL, ffOPTRD },
270 { efPDB, "-oq", "bfac", ffOPTWR },
271 { efPDB, "-ox", "xaver", ffOPTWR },
272 { efXVG, "-o", "rmsf", ffWRITE },
273 { efXVG, "-od", "rmsdev", ffOPTWR },
274 { efXVG, "-oc", "correl", ffOPTWR },
275 { efLOG, "-dir", "rmsf", ffOPTWR }
277 #define NFILE asize(fnm)
279 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW,
280 NFILE, fnm, asize(pargs), pargs, asize(desc), desc, 0, NULL,
286 bReadPDB = ftp2bSet(efPDB, NFILE, fnm);
287 devfn = opt2fn_null("-od", NFILE, fnm);
288 dirfn = opt2fn_null("-dir", NFILE, fnm);
290 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xref, NULL, box, TRUE);
291 snew(w_rls, top.atoms.nr);
293 fprintf(stderr, "Select group(s) for root mean square calculation\n");
294 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpnames);
297 for (i = 0; i < isize; i++)
299 w_rls[index[i]] = top.atoms.atom[index[i]].m;
302 /* Malloc the rmsf arrays */
303 snew(xav, isize*DIM);
305 for (i = 0; i < isize; i++)
317 get_stx_coordnum(opt2fn("-q", NFILE, fnm), &npdbatoms);
320 init_t_atoms(pdbatoms, npdbatoms, TRUE);
321 init_t_atoms(refatoms, npdbatoms, TRUE);
322 snew(pdbx, npdbatoms);
323 /* Read coordinates twice */
324 read_stx_conf(opt2fn("-q", NFILE, fnm), title, pdbatoms, pdbx, NULL, NULL, pdbbox);
325 read_stx_conf(opt2fn("-q", NFILE, fnm), title, refatoms, pdbx, NULL, NULL, pdbbox);
329 pdbatoms = &top.atoms;
330 refatoms = &top.atoms;
332 npdbatoms = pdbatoms->nr;
333 snew(pdbatoms->pdbinfo, npdbatoms);
334 copy_mat(box, pdbbox);
339 sub_xcm(xref, isize, index, top.atoms.atom, xcm, FALSE);
342 natom = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
346 gpbc = gmx_rmpbc_init(&top.idef, ePBC, natom);
349 /* Now read the trj again to compute fluctuations */
355 /* Remove periodic boundary */
356 gmx_rmpbc(gpbc, natom, box, x);
358 /* Set center of mass to zero */
359 sub_xcm(x, isize, index, top.atoms.atom, xcm, FALSE);
361 /* Fit to reference structure */
362 do_fit(natom, w_rls, xref, x);
365 /* Calculate Anisotropic U Tensor */
366 for (i = 0; i < isize; i++)
369 for (d = 0; d < DIM; d++)
371 xav[i*DIM + d] += x[aid][d];
372 for (m = 0; m < DIM; m++)
374 U[i][d*DIM + m] += x[aid][d]*x[aid][m];
381 /* Calculate RMS Deviation */
382 for (i = 0; (i < isize); i++)
385 for (d = 0; (d < DIM); d++)
387 rmsd_x[i][d] += sqr(x[aid][d]-xref[aid][d]);
394 while (read_next_x(oenv, status, &t, x, box));
399 gmx_rmpbc_done(gpbc);
403 invcount = 1.0/count;
404 snew(Uaver, DIM*DIM);
406 for (i = 0; i < isize; i++)
408 for (d = 0; d < DIM; d++)
410 xav[i*DIM + d] *= invcount;
412 for (d = 0; d < DIM; d++)
414 for (m = 0; m < DIM; m++)
416 U[i][d*DIM + m] = U[i][d*DIM + m]*invcount
417 - xav[i*DIM + d]*xav[i*DIM + m];
418 Uaver[3*d+m] += top.atoms.atom[index[i]].m*U[i][d*DIM + m];
421 totmass += top.atoms.atom[index[i]].m;
423 for (d = 0; d < DIM*DIM; d++)
430 for (d = 0; d < DIM*DIM; d++)
432 average_residues(NULL, U, d, isize, index, w_rls, &top.atoms);
438 for (i = 0; i < isize; i++)
441 pdbatoms->pdbinfo[aid].bAnisotropic = TRUE;
442 pdbatoms->pdbinfo[aid].uij[U11] = 1e6*U[i][XX*DIM + XX];
443 pdbatoms->pdbinfo[aid].uij[U22] = 1e6*U[i][YY*DIM + YY];
444 pdbatoms->pdbinfo[aid].uij[U33] = 1e6*U[i][ZZ*DIM + ZZ];
445 pdbatoms->pdbinfo[aid].uij[U12] = 1e6*U[i][XX*DIM + YY];
446 pdbatoms->pdbinfo[aid].uij[U13] = 1e6*U[i][XX*DIM + ZZ];
447 pdbatoms->pdbinfo[aid].uij[U23] = 1e6*U[i][YY*DIM + ZZ];
459 for (i = 0; i < isize; i++)
461 rmsf[i] = U[i][XX*DIM + XX] + U[i][YY*DIM + YY] + U[i][ZZ*DIM + ZZ];
466 fprintf(stdout, "\n");
467 print_dir(stdout, Uaver);
468 fp = gmx_ffopen(dirfn, "w");
469 print_dir(fp, Uaver);
473 for (i = 0; i < isize; i++)
479 /* Write RMSF output */
482 bfac = 8.0*M_PI*M_PI/3.0*100;
483 fp = xvgropen(ftp2fn(efXVG, NFILE, fnm), "B-Factors",
484 label, "(A\\b\\S\\So\\N\\S2\\N)", oenv);
485 xvgr_legend(fp, 2, leg, oenv);
486 for (i = 0; (i < isize); i++)
488 if (!bRes || i+1 == isize ||
489 top.atoms.atom[index[i]].resind != top.atoms.atom[index[i+1]].resind)
491 resind = top.atoms.atom[index[i]].resind;
492 pdb_bfac = find_pdb_bfac(pdbatoms, &top.atoms.resinfo[resind],
493 *(top.atoms.atomname[index[i]]));
495 fprintf(fp, "%5d %10.5f %10.5f\n",
496 bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : index[i]+1, rmsf[i]*bfac,
504 fp = xvgropen(ftp2fn(efXVG, NFILE, fnm), "RMS fluctuation", label, "(nm)", oenv);
505 for (i = 0; i < isize; i++)
507 if (!bRes || i+1 == isize ||
508 top.atoms.atom[index[i]].resind != top.atoms.atom[index[i+1]].resind)
510 fprintf(fp, "%5d %8.4f\n",
511 bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : index[i]+1, sqrt(rmsf[i]));
517 for (i = 0; i < isize; i++)
519 pdbatoms->pdbinfo[index[i]].bfac = 800*M_PI*M_PI/3.0*rmsf[i];
524 for (i = 0; i < isize; i++)
526 rmsf[i] = (rmsd_x[i][XX]+rmsd_x[i][YY]+rmsd_x[i][ZZ])/count;
530 average_residues(rmsf, NULL, 0, isize, index, w_rls, &top.atoms);
532 /* Write RMSD output */
533 fp = xvgropen(devfn, "RMS Deviation", label, "(nm)", oenv);
534 for (i = 0; i < isize; i++)
536 if (!bRes || i+1 == isize ||
537 top.atoms.atom[index[i]].resind != top.atoms.atom[index[i+1]].resind)
539 fprintf(fp, "%5d %8.4f\n",
540 bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : index[i]+1, sqrt(rmsf[i]));
546 if (opt2bSet("-oq", NFILE, fnm))
548 /* Write a .pdb file with B-factors and optionally anisou records */
549 for (i = 0; i < isize; i++)
551 rvec_inc(xref[index[i]], xcm);
553 write_sto_conf_indexed(opt2fn("-oq", NFILE, fnm), title, pdbatoms, pdbx,
554 NULL, ePBC, pdbbox, isize, index);
556 if (opt2bSet("-ox", NFILE, fnm))
558 /* Misuse xref as a temporary array */
559 for (i = 0; i < isize; i++)
561 for (d = 0; d < DIM; d++)
563 xref[index[i]][d] = xcm[d] + xav[i*DIM + d];
566 /* Write a .pdb file with B-factors and optionally anisou records */
567 write_sto_conf_indexed(opt2fn("-ox", NFILE, fnm), title, pdbatoms, xref, NULL,
568 ePBC, pdbbox, isize, index);
572 correlate_aniso(opt2fn("-oc", NFILE, fnm), refatoms, pdbatoms, oenv);
573 do_view(oenv, opt2fn("-oc", NFILE, fnm), "-nxy");
575 do_view(oenv, opt2fn("-o", NFILE, fnm), "-nxy");
578 do_view(oenv, opt2fn("-od", NFILE, fnm), "-nxy");