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,2016,2017,2018,2019, 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.
43 #include "gromacs/commandline/pargs.h"
44 #include "gromacs/commandline/viewit.h"
45 #include "gromacs/fileio/confio.h"
46 #include "gromacs/fileio/pdbio.h"
47 #include "gromacs/fileio/trxio.h"
48 #include "gromacs/fileio/xvgr.h"
49 #include "gromacs/gmxana/gmx_ana.h"
50 #include "gromacs/gmxana/princ.h"
51 #include "gromacs/linearalgebra/eigensolver.h"
52 #include "gromacs/math/do_fit.h"
53 #include "gromacs/math/functions.h"
54 #include "gromacs/math/utilities.h"
55 #include "gromacs/math/vec.h"
56 #include "gromacs/pbcutil/rmpbc.h"
57 #include "gromacs/topology/index.h"
58 #include "gromacs/topology/topology.h"
59 #include "gromacs/utility/arraysize.h"
60 #include "gromacs/utility/cstringutil.h"
61 #include "gromacs/utility/futil.h"
62 #include "gromacs/utility/smalloc.h"
64 static real find_pdb_bfac(const t_atoms* atoms, t_resinfo* ri, char* atomnm)
69 std::strcpy(rresnm, *ri->name);
71 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 && (std::strcmp(*atoms->resinfo[atoms->atom[i].resind].name, rresnm) == 0)
76 && (std::strstr(*atoms->atomname[i], atomnm) != nullptr))
83 fprintf(stderr, "\rCan not find %s%d-%s in pdbfile\n", rresnm, ri->nr, atomnm);
88 return atoms->pdbinfo[i].bfac;
91 static void correlate_aniso(const char* fn, t_atoms* ref, t_atoms* calc, const gmx_output_env_t* oenv)
96 fp = xvgropen(fn, "Correlation between X-Ray and Computed Uij", "X-Ray", "Computed", oenv);
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[],
116 const t_atoms* atoms)
128 for (i = 0; i < isize; i++)
130 av += w_rls[index[i]] * (f != nullptr ? f[i] : U[i][uind]);
131 m += w_rls[index[i]];
132 if (i + 1 == isize || atoms->atom[index[i]].resind != atoms->atom[index[i + 1]].resind)
137 for (j = start; j <= i; j++)
144 for (j = start; j <= i; j++)
156 static void print_dir(FILE* fp, real* Uaver)
158 real eigvec[DIM * DIM];
163 fprintf(fp, "MSF X Y Z\n");
164 for (d = 0; d < DIM; d++)
166 fprintf(fp, " %c ", 'X' + d - XX);
167 for (m = 0; m < DIM; m++)
169 fprintf(fp, " %9.2e", Uaver[3 * m + d]);
171 fprintf(fp, "%s\n", m == DIM ? " (nm^2)" : "");
174 for (m = 0; m < DIM * DIM; m++)
180 eigensolver(tmp, DIM, 0, DIM, eigval, eigvec);
182 fprintf(fp, "\n Eigenvectors\n\n");
183 fprintf(fp, "Eigv %-8.2e %-8.2e %-8.2e (nm^2)\n\n", eigval[2], eigval[1], eigval[0]);
184 for (d = 0; d < DIM; d++)
186 fprintf(fp, " %c ", 'X' + d - XX);
187 for (m = DIM - 1; m >= 0; m--)
189 fprintf(fp, "%7.4f ", eigvec[3 * m + d]);
195 int gmx_rmsf(int argc, char* argv[])
197 const char* desc[] = {
198 "[THISMODULE] computes the root mean square fluctuation (RMSF, i.e. standard ",
199 "deviation) of atomic positions in the trajectory (supplied with [TT]-f[tt])",
200 "after (optionally) fitting to a reference frame (supplied with [TT]-s[tt]).[PAR]",
201 "With option [TT]-oq[tt] the RMSF values are converted to B-factor",
202 "values, which are written to a [REF].pdb[ref] file. By default, the coordinates",
203 "in this output file are taken from the structure file provided with [TT]-s[tt],"
204 "although you can also use coordinates read from a different [REF].pdb[ref] file"
205 "provided with [TT]-q[tt]. There is very little error checking, so in this case"
206 "it is your responsibility to make sure all atoms in the structure file"
207 "and [REF].pdb[ref] file correspond exactly to each other.[PAR]",
208 "Option [TT]-ox[tt] writes the B-factors to a file with the average",
209 "coordinates in the trajectory.[PAR]",
210 "With the option [TT]-od[tt] the root mean square deviation with",
211 "respect to the reference structure is calculated.[PAR]",
212 "With the option [TT]-aniso[tt], [THISMODULE] will compute anisotropic",
213 "temperature factors and then it will also output average coordinates",
214 "and a [REF].pdb[ref] file with ANISOU records (corresonding to the [TT]-oq[tt]",
215 "or [TT]-ox[tt] option). Please note that the U values",
216 "are orientation-dependent, so before comparison with experimental data",
217 "you should verify that you fit to the experimental coordinates.[PAR]",
218 "When a [REF].pdb[ref] input file is passed to the program and the [TT]-aniso[tt]",
220 "a correlation plot of the Uij will be created, if any anisotropic",
221 "temperature factors are present in the [REF].pdb[ref] file.[PAR]",
222 "With option [TT]-dir[tt] the average MSF (3x3) matrix is diagonalized.",
223 "This shows the directions in which the atoms fluctuate the most and",
226 static gmx_bool bRes = FALSE, bAniso = FALSE, bFit = TRUE;
228 { "-res", FALSE, etBOOL, { &bRes }, "Calculate averages for each residue" },
229 { "-aniso", FALSE, etBOOL, { &bAniso }, "Compute anisotropic termperature factors" },
234 "Do a least squares superposition before computing RMSF. Without this you must "
235 "make sure that the reference structure and the trajectory match." }
238 int i, m, teller = 0;
243 t_atoms * pdbatoms, *refatoms;
246 rvec * x, *pdbx, *xref;
250 FILE* fp; /* the graphics file */
251 const char *devfn, *dirfn;
259 real bfac, pdb_bfac, *Uaver;
262 rvec* rmsd_x = nullptr;
263 double * rmsf, invcount, totmass;
267 gmx_rmpbc_t gpbc = nullptr;
269 gmx_output_env_t* oenv;
271 const char* leg[2] = { "MD", "X-Ray" };
273 t_filenm fnm[] = { { efTRX, "-f", nullptr, ffREAD }, { efTPS, nullptr, nullptr, ffREAD },
274 { efNDX, nullptr, nullptr, ffOPTRD }, { efPDB, "-q", nullptr, ffOPTRD },
275 { efPDB, "-oq", "bfac", ffOPTWR }, { efPDB, "-ox", "xaver", ffOPTWR },
276 { efXVG, "-o", "rmsf", ffWRITE }, { efXVG, "-od", "rmsdev", ffOPTWR },
277 { efXVG, "-oc", "correl", ffOPTWR }, { efLOG, "-dir", "rmsf", ffOPTWR } };
278 #define NFILE asize(fnm)
280 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW, NFILE, fnm, asize(pargs),
281 pargs, asize(desc), desc, 0, nullptr, &oenv))
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), &top, &ePBC, &xref, nullptr, box, TRUE);
291 const char* title = *top.name;
292 snew(w_rls, top.atoms.nr);
294 fprintf(stderr, "Select group(s) for root mean square calculation\n");
295 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpnames);
298 for (i = 0; i < isize; i++)
300 w_rls[index[i]] = top.atoms.atom[index[i]].m;
303 /* Malloc the rmsf arrays */
304 snew(xav, isize * DIM);
306 for (i = 0; i < isize; i++)
308 snew(U[i], DIM * DIM);
320 /* Read coordinates twice */
321 read_tps_conf(opt2fn("-q", NFILE, fnm), top_pdb, nullptr, nullptr, nullptr, pdbbox, FALSE);
323 *pdbatoms = top_pdb->atoms;
324 read_tps_conf(opt2fn("-q", NFILE, fnm), top_pdb, nullptr, &pdbx, nullptr, pdbbox, FALSE);
325 /* TODO Should this assert that top_pdb->atoms.nr == top.atoms.nr?
326 * See discussion at https://gerrit.gromacs.org/#/c/6430/1 */
327 title = *top_pdb->name;
329 *refatoms = top_pdb->atoms;
334 pdbatoms = &top.atoms;
335 refatoms = &top.atoms;
337 snew(pdbatoms->pdbinfo, pdbatoms->nr);
338 pdbatoms->havePdbInfo = TRUE;
339 copy_mat(box, pdbbox);
344 sub_xcm(xref, isize, index, top.atoms.atom, xcm, FALSE);
347 natom = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
351 gpbc = gmx_rmpbc_init(&top.idef, ePBC, natom);
354 /* Now read the trj again to compute fluctuations */
360 /* Remove periodic boundary */
361 gmx_rmpbc(gpbc, natom, box, x);
363 /* Set center of mass to zero */
364 sub_xcm(x, isize, index, top.atoms.atom, xcm, FALSE);
366 /* Fit to reference structure */
367 do_fit(natom, w_rls, xref, x);
370 /* Calculate Anisotropic U Tensor */
371 for (i = 0; i < isize; i++)
374 for (d = 0; d < DIM; d++)
376 xav[i * DIM + d] += x[aid][d];
377 for (m = 0; m < DIM; m++)
379 U[i][d * DIM + m] += x[aid][d] * x[aid][m];
386 /* Calculate RMS Deviation */
387 for (i = 0; (i < isize); i++)
390 for (d = 0; (d < DIM); d++)
392 rmsd_x[i][d] += gmx::square(x[aid][d] - xref[aid][d]);
398 } while (read_next_x(oenv, status, &t, x, box));
403 gmx_rmpbc_done(gpbc);
407 invcount = 1.0 / count;
408 snew(Uaver, DIM * DIM);
410 for (i = 0; i < isize; i++)
412 for (d = 0; d < DIM; d++)
414 xav[i * DIM + d] *= invcount;
416 for (d = 0; d < DIM; d++)
418 for (m = 0; m < DIM; m++)
420 U[i][d * DIM + m] = U[i][d * DIM + m] * invcount - xav[i * DIM + d] * xav[i * DIM + m];
421 Uaver[3 * d + m] += top.atoms.atom[index[i]].m * U[i][d * DIM + m];
424 totmass += top.atoms.atom[index[i]].m;
426 for (d = 0; d < DIM * DIM; d++)
433 for (d = 0; d < DIM * DIM; d++)
435 average_residues(nullptr, U, d, isize, index, w_rls, &top.atoms);
441 for (i = 0; i < isize; i++)
444 pdbatoms->pdbinfo[aid].bAnisotropic = TRUE;
445 pdbatoms->pdbinfo[aid].uij[U11] = static_cast<int>(1e6 * U[i][XX * DIM + XX]);
446 pdbatoms->pdbinfo[aid].uij[U22] = static_cast<int>(1e6 * U[i][YY * DIM + YY]);
447 pdbatoms->pdbinfo[aid].uij[U33] = static_cast<int>(1e6 * U[i][ZZ * DIM + ZZ]);
448 pdbatoms->pdbinfo[aid].uij[U12] = static_cast<int>(1e6 * U[i][XX * DIM + YY]);
449 pdbatoms->pdbinfo[aid].uij[U13] = static_cast<int>(1e6 * U[i][XX * DIM + ZZ]);
450 pdbatoms->pdbinfo[aid].uij[U23] = static_cast<int>(1e6 * U[i][YY * DIM + ZZ]);
462 for (i = 0; i < isize; i++)
464 rmsf[i] = U[i][XX * DIM + XX] + U[i][YY * DIM + YY] + U[i][ZZ * DIM + ZZ];
469 fprintf(stdout, "\n");
470 print_dir(stdout, Uaver);
471 fp = gmx_ffopen(dirfn, "w");
472 print_dir(fp, Uaver);
476 for (i = 0; i < isize; i++)
482 /* Write RMSF output */
485 bfac = 8.0 * M_PI * M_PI / 3.0 * 100;
486 fp = xvgropen(ftp2fn(efXVG, NFILE, fnm), "B-Factors", label, "(A\\b\\S\\So\\N\\S2\\N)", oenv);
487 xvgr_legend(fp, 2, leg, oenv);
488 for (i = 0; (i < isize); i++)
490 if (!bRes || i + 1 == isize
491 || top.atoms.atom[index[i]].resind != top.atoms.atom[index[i + 1]].resind)
493 resind = top.atoms.atom[index[i]].resind;
494 pdb_bfac = find_pdb_bfac(pdbatoms, &top.atoms.resinfo[resind],
495 *(top.atoms.atomname[index[i]]));
497 fprintf(fp, "%5d %10.5f %10.5f\n",
498 bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : index[i] + 1,
499 rmsf[i] * bfac, pdb_bfac);
506 fp = xvgropen(ftp2fn(efXVG, NFILE, fnm), "RMS fluctuation", label, "(nm)", oenv);
507 for (i = 0; i < isize; i++)
509 if (!bRes || i + 1 == isize
510 || top.atoms.atom[index[i]].resind != top.atoms.atom[index[i + 1]].resind)
512 fprintf(fp, "%5d %8.4f\n",
513 bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : index[i] + 1,
520 for (i = 0; i < isize; i++)
522 pdbatoms->pdbinfo[index[i]].bfac = 800 * M_PI * M_PI / 3.0 * rmsf[i];
527 for (i = 0; i < isize; i++)
529 rmsf[i] = (rmsd_x[i][XX] + rmsd_x[i][YY] + rmsd_x[i][ZZ]) / count;
533 average_residues(rmsf, nullptr, 0, isize, index, w_rls, &top.atoms);
535 /* Write RMSD output */
536 fp = xvgropen(devfn, "RMS Deviation", label, "(nm)", oenv);
537 for (i = 0; i < isize; i++)
539 if (!bRes || i + 1 == isize
540 || top.atoms.atom[index[i]].resind != top.atoms.atom[index[i + 1]].resind)
542 fprintf(fp, "%5d %8.4f\n",
543 bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : index[i] + 1,
550 if (opt2bSet("-oq", NFILE, fnm))
552 /* Write a .pdb file with B-factors and optionally anisou records */
553 for (i = 0; i < isize; i++)
555 rvec_inc(pdbx[index[i]], xcm);
557 write_sto_conf_indexed(opt2fn("-oq", NFILE, fnm), title, pdbatoms, pdbx, nullptr, ePBC,
558 pdbbox, isize, index);
560 if (opt2bSet("-ox", NFILE, fnm))
563 snew(bFactorX, top.atoms.nr);
564 for (i = 0; i < isize; i++)
566 for (d = 0; d < DIM; d++)
568 bFactorX[index[i]][d] = xcm[d] + xav[i * DIM + d];
571 /* Write a .pdb file with B-factors and optionally anisou records */
572 write_sto_conf_indexed(opt2fn("-ox", NFILE, fnm), title, pdbatoms, bFactorX, nullptr, ePBC,
573 pdbbox, isize, index);
578 correlate_aniso(opt2fn("-oc", NFILE, fnm), refatoms, pdbatoms, oenv);
579 do_view(oenv, opt2fn("-oc", NFILE, fnm), "-nxy");
581 do_view(oenv, opt2fn("-o", NFILE, fnm), "-nxy");
584 do_view(oenv, opt2fn("-od", NFILE, fnm), "-nxy");