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,2015, 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.
42 #include "gromacs/commandline/pargs.h"
43 #include "gromacs/fileio/confio.h"
44 #include "gromacs/fileio/pdbio.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/cstringutil.h"
58 #include "gromacs/utility/futil.h"
59 #include "gromacs/utility/smalloc.h"
61 static real find_pdb_bfac(t_atoms *atoms, t_resinfo *ri, char *atomnm)
66 std::strcpy(rresnm, *ri->name);
68 for (i = 0; (i < atoms->nr); i++)
70 if ((ri->nr == atoms->resinfo[atoms->atom[i].resind].nr) &&
71 (ri->ic == atoms->resinfo[atoms->atom[i].resind].ic) &&
72 (std::strcmp(*atoms->resinfo[atoms->atom[i].resind].name, rresnm) == 0) &&
73 (std::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++)
98 if (ref->pdbinfo[i].bAnisotropic)
100 for (j = U11; (j <= U23); j++)
102 fprintf(fp, "%10d %10d\n", ref->pdbinfo[i].uij[j], calc->pdbinfo[i].uij[j]);
109 static void average_residues(double f[], double **U, int uind,
110 int isize, atom_id index[], real w_rls[],
119 for (i = 0; i < isize; i++)
121 av += w_rls[index[i]]*(f != NULL ? f[i] : U[i][uind]);
122 m += w_rls[index[i]];
124 atoms->atom[index[i]].resind != atoms->atom[index[i+1]].resind)
129 for (j = start; j <= i; j++)
136 for (j = start; j <= i; j++)
148 void print_dir(FILE *fp, real *Uaver)
150 real eigvec[DIM*DIM];
155 fprintf(fp, "MSF X Y Z\n");
156 for (d = 0; d < DIM; d++)
158 fprintf(fp, " %c ", 'X'+d-XX);
159 for (m = 0; m < DIM; m++)
161 fprintf(fp, " %9.2e", Uaver[3*m+d]);
163 fprintf(fp, "%s\n", m == DIM ? " (nm^2)" : "");
166 for (m = 0; m < DIM*DIM; m++)
172 eigensolver(tmp, DIM, 0, DIM, eigval, eigvec);
174 fprintf(fp, "\n Eigenvectors\n\n");
175 fprintf(fp, "Eigv %-8.2e %-8.2e %-8.2e (nm^2)\n\n",
176 eigval[2], eigval[1], eigval[0]);
177 for (d = 0; d < DIM; d++)
179 fprintf(fp, " %c ", 'X'+d-XX);
180 for (m = DIM-1; m >= 0; m--)
182 fprintf(fp, "%7.4f ", eigvec[3*m+d]);
188 int gmx_rmsf(int argc, char *argv[])
190 const char *desc[] = {
191 "[THISMODULE] computes the root mean square fluctuation (RMSF, i.e. standard ",
192 "deviation) of atomic positions in the trajectory (supplied with [TT]-f[tt])",
193 "after (optionally) fitting to a reference frame (supplied with [TT]-s[tt]).[PAR]",
194 "With option [TT]-oq[tt] the RMSF values are converted to B-factor",
195 "values, which are written to a [REF].pdb[ref] file with the coordinates, of the",
196 "structure file, or of a [REF].pdb[ref] file when [TT]-q[tt] is specified.",
197 "Option [TT]-ox[tt] writes the B-factors to a file with the average",
199 "With the option [TT]-od[tt] the root mean square deviation with",
200 "respect to the reference structure is calculated.[PAR]",
201 "With the option [TT]-aniso[tt], [THISMODULE] will compute anisotropic",
202 "temperature factors and then it will also output average coordinates",
203 "and a [REF].pdb[ref] file with ANISOU records (corresonding to the [TT]-oq[tt]",
204 "or [TT]-ox[tt] option). Please note that the U values",
205 "are orientation-dependent, so before comparison with experimental data",
206 "you should verify that you fit to the experimental coordinates.[PAR]",
207 "When a [REF].pdb[ref] input file is passed to the program and the [TT]-aniso[tt]",
209 "a correlation plot of the Uij will be created, if any anisotropic",
210 "temperature factors are present in the [REF].pdb[ref] file.[PAR]",
211 "With option [TT]-dir[tt] the average MSF (3x3) matrix is diagonalized.",
212 "This shows the directions in which the atoms fluctuate the most and",
215 static gmx_bool bRes = FALSE, bAniso = FALSE, bFit = TRUE;
217 { "-res", FALSE, etBOOL, {&bRes},
218 "Calculate averages for each residue" },
219 { "-aniso", FALSE, etBOOL, {&bAniso},
220 "Compute anisotropic termperature factors" },
221 { "-fit", FALSE, etBOOL, {&bFit},
222 "Do a least squares superposition before computing RMSF. Without this you must make sure that the reference structure and the trajectory match." }
225 int i, m, teller = 0;
230 t_atoms *pdbatoms, *refatoms;
233 rvec *x, *pdbx, *xref;
239 FILE *fp; /* the graphics file */
240 const char *devfn, *dirfn;
248 real bfac, pdb_bfac, *Uaver;
252 double *rmsf, invcount, totmass;
256 gmx_rmpbc_t gpbc = NULL;
260 const char *leg[2] = { "MD", "X-Ray" };
263 { efTRX, "-f", NULL, ffREAD },
264 { efTPS, NULL, NULL, ffREAD },
265 { efNDX, NULL, NULL, ffOPTRD },
266 { efPDB, "-q", NULL, ffOPTRD },
267 { efPDB, "-oq", "bfac", ffOPTWR },
268 { efPDB, "-ox", "xaver", ffOPTWR },
269 { efXVG, "-o", "rmsf", ffWRITE },
270 { efXVG, "-od", "rmsdev", ffOPTWR },
271 { efXVG, "-oc", "correl", ffOPTWR },
272 { efLOG, "-dir", "rmsf", ffOPTWR }
274 #define NFILE asize(fnm)
276 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW,
277 NFILE, fnm, asize(pargs), pargs, asize(desc), desc, 0, NULL,
283 bReadPDB = ftp2bSet(efPDB, NFILE, fnm);
284 devfn = opt2fn_null("-od", NFILE, fnm);
285 dirfn = opt2fn_null("-dir", NFILE, fnm);
287 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xref, NULL, box, TRUE);
288 snew(w_rls, top.atoms.nr);
290 fprintf(stderr, "Select group(s) for root mean square calculation\n");
291 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpnames);
294 for (i = 0; i < isize; i++)
296 w_rls[index[i]] = top.atoms.atom[index[i]].m;
299 /* Malloc the rmsf arrays */
300 snew(xav, isize*DIM);
302 for (i = 0; i < isize; i++)
314 get_stx_coordnum(opt2fn("-q", NFILE, fnm), &npdbatoms);
317 init_t_atoms(pdbatoms, npdbatoms, TRUE);
318 init_t_atoms(refatoms, npdbatoms, TRUE);
319 snew(pdbx, npdbatoms);
320 /* Read coordinates twice */
321 read_stx_conf(opt2fn("-q", NFILE, fnm), title, pdbatoms, pdbx, NULL, NULL, pdbbox);
322 read_stx_conf(opt2fn("-q", NFILE, fnm), title, refatoms, pdbx, NULL, NULL, pdbbox);
326 pdbatoms = &top.atoms;
327 refatoms = &top.atoms;
329 npdbatoms = pdbatoms->nr;
330 snew(pdbatoms->pdbinfo, npdbatoms);
331 copy_mat(box, pdbbox);
336 sub_xcm(xref, isize, index, top.atoms.atom, xcm, FALSE);
339 natom = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
343 gpbc = gmx_rmpbc_init(&top.idef, ePBC, natom);
346 /* Now read the trj again to compute fluctuations */
352 /* Remove periodic boundary */
353 gmx_rmpbc(gpbc, natom, box, x);
355 /* Set center of mass to zero */
356 sub_xcm(x, isize, index, top.atoms.atom, xcm, FALSE);
358 /* Fit to reference structure */
359 do_fit(natom, w_rls, xref, x);
362 /* Calculate Anisotropic U Tensor */
363 for (i = 0; i < isize; i++)
366 for (d = 0; d < DIM; d++)
368 xav[i*DIM + d] += x[aid][d];
369 for (m = 0; m < DIM; m++)
371 U[i][d*DIM + m] += x[aid][d]*x[aid][m];
378 /* Calculate RMS Deviation */
379 for (i = 0; (i < isize); i++)
382 for (d = 0; (d < DIM); d++)
384 rmsd_x[i][d] += sqr(x[aid][d]-xref[aid][d]);
391 while (read_next_x(oenv, status, &t, x, box));
396 gmx_rmpbc_done(gpbc);
400 invcount = 1.0/count;
401 snew(Uaver, DIM*DIM);
403 for (i = 0; i < isize; i++)
405 for (d = 0; d < DIM; d++)
407 xav[i*DIM + d] *= invcount;
409 for (d = 0; d < DIM; d++)
411 for (m = 0; m < DIM; m++)
413 U[i][d*DIM + m] = U[i][d*DIM + m]*invcount
414 - xav[i*DIM + d]*xav[i*DIM + m];
415 Uaver[3*d+m] += top.atoms.atom[index[i]].m*U[i][d*DIM + m];
418 totmass += top.atoms.atom[index[i]].m;
420 for (d = 0; d < DIM*DIM; d++)
427 for (d = 0; d < DIM*DIM; d++)
429 average_residues(NULL, U, d, isize, index, w_rls, &top.atoms);
435 for (i = 0; i < isize; i++)
438 pdbatoms->pdbinfo[aid].bAnisotropic = TRUE;
439 pdbatoms->pdbinfo[aid].uij[U11] = static_cast<int>(1e6*U[i][XX*DIM + XX]);
440 pdbatoms->pdbinfo[aid].uij[U22] = static_cast<int>(1e6*U[i][YY*DIM + YY]);
441 pdbatoms->pdbinfo[aid].uij[U33] = static_cast<int>(1e6*U[i][ZZ*DIM + ZZ]);
442 pdbatoms->pdbinfo[aid].uij[U12] = static_cast<int>(1e6*U[i][XX*DIM + YY]);
443 pdbatoms->pdbinfo[aid].uij[U13] = static_cast<int>(1e6*U[i][XX*DIM + ZZ]);
444 pdbatoms->pdbinfo[aid].uij[U23] = static_cast<int>(1e6*U[i][YY*DIM + ZZ]);
456 for (i = 0; i < isize; i++)
458 rmsf[i] = U[i][XX*DIM + XX] + U[i][YY*DIM + YY] + U[i][ZZ*DIM + ZZ];
463 fprintf(stdout, "\n");
464 print_dir(stdout, Uaver);
465 fp = gmx_ffopen(dirfn, "w");
466 print_dir(fp, Uaver);
470 for (i = 0; i < isize; i++)
476 /* Write RMSF output */
479 bfac = 8.0*M_PI*M_PI/3.0*100;
480 fp = xvgropen(ftp2fn(efXVG, NFILE, fnm), "B-Factors",
481 label, "(A\\b\\S\\So\\N\\S2\\N)", oenv);
482 xvgr_legend(fp, 2, leg, oenv);
483 for (i = 0; (i < isize); i++)
485 if (!bRes || i+1 == isize ||
486 top.atoms.atom[index[i]].resind != top.atoms.atom[index[i+1]].resind)
488 resind = top.atoms.atom[index[i]].resind;
489 pdb_bfac = find_pdb_bfac(pdbatoms, &top.atoms.resinfo[resind],
490 *(top.atoms.atomname[index[i]]));
492 fprintf(fp, "%5d %10.5f %10.5f\n",
493 bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : index[i]+1, rmsf[i]*bfac,
501 fp = xvgropen(ftp2fn(efXVG, NFILE, fnm), "RMS fluctuation", label, "(nm)", oenv);
502 for (i = 0; i < isize; i++)
504 if (!bRes || i+1 == isize ||
505 top.atoms.atom[index[i]].resind != top.atoms.atom[index[i+1]].resind)
507 fprintf(fp, "%5d %8.4f\n",
508 bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : index[i]+1, std::sqrt(rmsf[i]));
514 for (i = 0; i < isize; i++)
516 pdbatoms->pdbinfo[index[i]].bfac = 800*M_PI*M_PI/3.0*rmsf[i];
521 for (i = 0; i < isize; i++)
523 rmsf[i] = (rmsd_x[i][XX]+rmsd_x[i][YY]+rmsd_x[i][ZZ])/count;
527 average_residues(rmsf, NULL, 0, isize, index, w_rls, &top.atoms);
529 /* Write RMSD output */
530 fp = xvgropen(devfn, "RMS Deviation", label, "(nm)", oenv);
531 for (i = 0; i < isize; i++)
533 if (!bRes || i+1 == isize ||
534 top.atoms.atom[index[i]].resind != top.atoms.atom[index[i+1]].resind)
536 fprintf(fp, "%5d %8.4f\n",
537 bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : index[i]+1, std::sqrt(rmsf[i]));
543 if (opt2bSet("-oq", NFILE, fnm))
545 /* Write a .pdb file with B-factors and optionally anisou records */
546 for (i = 0; i < isize; i++)
548 rvec_inc(xref[index[i]], xcm);
550 write_sto_conf_indexed(opt2fn("-oq", NFILE, fnm), title, pdbatoms, pdbx,
551 NULL, ePBC, pdbbox, isize, index);
553 if (opt2bSet("-ox", NFILE, fnm))
555 /* Misuse xref as a temporary array */
556 for (i = 0; i < isize; i++)
558 for (d = 0; d < DIM; d++)
560 xref[index[i]][d] = xcm[d] + xav[i*DIM + d];
563 /* Write a .pdb file with B-factors and optionally anisou records */
564 write_sto_conf_indexed(opt2fn("-ox", NFILE, fnm), title, pdbatoms, xref, NULL,
565 ePBC, pdbbox, isize, index);
569 correlate_aniso(opt2fn("-oc", NFILE, fnm), refatoms, pdbatoms, oenv);
570 do_view(oenv, opt2fn("-oc", NFILE, fnm), "-nxy");
572 do_view(oenv, opt2fn("-o", NFILE, fnm), "-nxy");
575 do_view(oenv, opt2fn("-od", NFILE, fnm), "-nxy");