Tidy: modernize-use-nullptr
[gromacs.git] / src / gromacs / gmxana / gmx_rmsf.cpp
blob9c3989a279ecb0018b11368faa684bb3e9a8234d
1 /*
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, 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.
37 #include "gmxpre.h"
39 #include <cassert>
40 #include <cmath>
41 #include <cstring>
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)
66 char rresnm[8];
67 int i;
69 std::strcpy(rresnm, *ri->name);
70 rresnm[3] = '\0';
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))
78 break;
81 if (i == atoms->nr)
83 fprintf(stderr, "\rCan not find %s%d-%s in pdbfile\n",
84 rresnm, ri->nr, atomnm);
85 fflush(stderr);
86 return 0.0;
89 return atoms->pdbinfo[i].bfac;
92 void correlate_aniso(const char *fn, t_atoms *ref, t_atoms *calc,
93 const gmx_output_env_t *oenv)
95 FILE *fp;
96 int i, j;
98 fp = xvgropen(fn, "Correlation between X-Ray and Computed Uij", "X-Ray",
99 "Computed", oenv);
100 for (i = 0; (i < ref->nr); i++)
102 if (ref->pdbinfo[i].bAnisotropic)
104 for (j = U11; (j <= U23); j++)
106 fprintf(fp, "%10d %10d\n", ref->pdbinfo[i].uij[j], calc->pdbinfo[i].uij[j]);
110 xvgrclose(fp);
113 static void average_residues(double f[], double **U, int uind,
114 int isize, int index[], real w_rls[],
115 const t_atoms *atoms)
117 int i, j, start;
118 double av, m;
120 start = 0;
121 av = 0;
122 m = 0;
123 if (!f)
125 assert(U);
127 for (i = 0; i < isize; i++)
129 av += w_rls[index[i]]*(f != nullptr ? f[i] : U[i][uind]);
130 m += w_rls[index[i]];
131 if (i+1 == isize ||
132 atoms->atom[index[i]].resind != atoms->atom[index[i+1]].resind)
134 av /= m;
135 if (f != nullptr)
137 for (j = start; j <= i; j++)
139 f[i] = av;
142 else
144 for (j = start; j <= i; j++)
146 U[j][uind] = av;
149 start = i+1;
150 av = 0;
151 m = 0;
156 void print_dir(FILE *fp, real *Uaver)
158 real eigvec[DIM*DIM];
159 real tmp[DIM*DIM];
160 rvec eigval;
161 int d, m;
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++)
176 tmp[m] = Uaver[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",
184 eigval[2], eigval[1], eigval[0]);
185 for (d = 0; d < DIM; d++)
187 fprintf(fp, " %c ", 'X'+d-XX);
188 for (m = DIM-1; m >= 0; m--)
190 fprintf(fp, "%7.4f ", eigvec[3*m+d]);
192 fprintf(fp, "\n");
196 int gmx_rmsf(int argc, char *argv[])
198 const char *desc[] = {
199 "[THISMODULE] computes the root mean square fluctuation (RMSF, i.e. standard ",
200 "deviation) of atomic positions in the trajectory (supplied with [TT]-f[tt])",
201 "after (optionally) fitting to a reference frame (supplied with [TT]-s[tt]).[PAR]",
202 "With option [TT]-oq[tt] the RMSF values are converted to B-factor",
203 "values, which are written to a [REF].pdb[ref] file with the coordinates, of the",
204 "structure file, or of a [REF].pdb[ref] file when [TT]-q[tt] is specified.",
205 "Option [TT]-ox[tt] writes the B-factors to a file with the average",
206 "coordinates.[PAR]",
207 "With the option [TT]-od[tt] the root mean square deviation with",
208 "respect to the reference structure is calculated.[PAR]",
209 "With the option [TT]-aniso[tt], [THISMODULE] will compute anisotropic",
210 "temperature factors and then it will also output average coordinates",
211 "and a [REF].pdb[ref] file with ANISOU records (corresonding to the [TT]-oq[tt]",
212 "or [TT]-ox[tt] option). Please note that the U values",
213 "are orientation-dependent, so before comparison with experimental data",
214 "you should verify that you fit to the experimental coordinates.[PAR]",
215 "When a [REF].pdb[ref] input file is passed to the program and the [TT]-aniso[tt]",
216 "flag is set",
217 "a correlation plot of the Uij will be created, if any anisotropic",
218 "temperature factors are present in the [REF].pdb[ref] file.[PAR]",
219 "With option [TT]-dir[tt] the average MSF (3x3) matrix is diagonalized.",
220 "This shows the directions in which the atoms fluctuate the most and",
221 "the least."
223 static gmx_bool bRes = FALSE, bAniso = FALSE, bFit = TRUE;
224 t_pargs pargs[] = {
225 { "-res", FALSE, etBOOL, {&bRes},
226 "Calculate averages for each residue" },
227 { "-aniso", FALSE, etBOOL, {&bAniso},
228 "Compute anisotropic termperature factors" },
229 { "-fit", FALSE, etBOOL, {&bFit},
230 "Do a least squares superposition before computing RMSF. Without this you must make sure that the reference structure and the trajectory match." }
232 int natom;
233 int i, m, teller = 0;
234 real t, *w_rls;
236 t_topology top;
237 int ePBC;
238 t_atoms *pdbatoms, *refatoms;
240 matrix box, pdbbox;
241 rvec *x, *pdbx, *xref;
242 t_trxstatus *status;
243 const char *label;
245 FILE *fp; /* the graphics file */
246 const char *devfn, *dirfn;
247 int resind;
249 gmx_bool bReadPDB;
250 int *index;
251 int isize;
252 char *grpnames;
254 real bfac, pdb_bfac, *Uaver;
255 double **U, *xav;
256 int aid;
257 rvec *rmsd_x = nullptr;
258 double *rmsf, invcount, totmass;
259 int d;
260 real count = 0;
261 rvec xcm;
262 gmx_rmpbc_t gpbc = nullptr;
264 gmx_output_env_t *oenv;
266 const char *leg[2] = { "MD", "X-Ray" };
268 t_filenm fnm[] = {
269 { efTRX, "-f", nullptr, ffREAD },
270 { efTPS, nullptr, nullptr, ffREAD },
271 { efNDX, nullptr, nullptr, ffOPTRD },
272 { efPDB, "-q", nullptr, ffOPTRD },
273 { efPDB, "-oq", "bfac", ffOPTWR },
274 { efPDB, "-ox", "xaver", ffOPTWR },
275 { efXVG, "-o", "rmsf", ffWRITE },
276 { efXVG, "-od", "rmsdev", ffOPTWR },
277 { efXVG, "-oc", "correl", ffOPTWR },
278 { efLOG, "-dir", "rmsf", ffOPTWR }
280 #define NFILE asize(fnm)
282 if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW,
283 NFILE, fnm, asize(pargs), pargs, asize(desc), desc, 0, nullptr,
284 &oenv))
286 return 0;
289 bReadPDB = ftp2bSet(efPDB, NFILE, fnm);
290 devfn = opt2fn_null("-od", NFILE, fnm);
291 dirfn = opt2fn_null("-dir", NFILE, fnm);
293 read_tps_conf(ftp2fn(efTPS, NFILE, fnm), &top, &ePBC, &xref, nullptr, box, TRUE);
294 const char *title = *top.name;
295 snew(w_rls, top.atoms.nr);
297 fprintf(stderr, "Select group(s) for root mean square calculation\n");
298 get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpnames);
300 /* Set the weight */
301 for (i = 0; i < isize; i++)
303 w_rls[index[i]] = top.atoms.atom[index[i]].m;
306 /* Malloc the rmsf arrays */
307 snew(xav, isize*DIM);
308 snew(U, isize);
309 for (i = 0; i < isize; i++)
311 snew(U[i], DIM*DIM);
313 snew(rmsf, isize);
314 if (devfn)
316 snew(rmsd_x, isize);
319 if (bReadPDB)
321 t_topology *top_pdb;
322 snew(top_pdb, 1);
323 /* Read coordinates twice */
324 read_tps_conf(opt2fn("-q", NFILE, fnm), top_pdb, nullptr, nullptr, nullptr, pdbbox, FALSE);
325 snew(pdbatoms, 1);
326 *pdbatoms = top_pdb->atoms;
327 read_tps_conf(opt2fn("-q", NFILE, fnm), top_pdb, nullptr, &pdbx, nullptr, pdbbox, FALSE);
328 title = *top_pdb->name;
329 snew(refatoms, 1);
330 *refatoms = top_pdb->atoms;
331 sfree(top_pdb);
333 else
335 pdbatoms = &top.atoms;
336 refatoms = &top.atoms;
337 pdbx = xref;
338 snew(pdbatoms->pdbinfo, pdbatoms->nr);
339 pdbatoms->havePdbInfo = TRUE;
340 copy_mat(box, pdbbox);
343 if (bFit)
345 sub_xcm(xref, isize, index, top.atoms.atom, xcm, FALSE);
348 natom = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
350 if (bFit)
352 gpbc = gmx_rmpbc_init(&top.idef, ePBC, natom);
355 /* Now read the trj again to compute fluctuations */
356 teller = 0;
359 if (bFit)
361 /* Remove periodic boundary */
362 gmx_rmpbc(gpbc, natom, box, x);
364 /* Set center of mass to zero */
365 sub_xcm(x, isize, index, top.atoms.atom, xcm, FALSE);
367 /* Fit to reference structure */
368 do_fit(natom, w_rls, xref, x);
371 /* Calculate Anisotropic U Tensor */
372 for (i = 0; i < isize; i++)
374 aid = index[i];
375 for (d = 0; d < DIM; d++)
377 xav[i*DIM + d] += x[aid][d];
378 for (m = 0; m < DIM; m++)
380 U[i][d*DIM + m] += x[aid][d]*x[aid][m];
385 if (devfn)
387 /* Calculate RMS Deviation */
388 for (i = 0; (i < isize); i++)
390 aid = index[i];
391 for (d = 0; (d < DIM); d++)
393 rmsd_x[i][d] += gmx::square(x[aid][d]-xref[aid][d]);
397 count += 1.0;
398 teller++;
400 while (read_next_x(oenv, status, &t, x, box));
401 close_trj(status);
403 if (bFit)
405 gmx_rmpbc_done(gpbc);
409 invcount = 1.0/count;
410 snew(Uaver, DIM*DIM);
411 totmass = 0;
412 for (i = 0; i < isize; i++)
414 for (d = 0; d < DIM; d++)
416 xav[i*DIM + d] *= invcount;
418 for (d = 0; d < DIM; d++)
420 for (m = 0; m < DIM; m++)
422 U[i][d*DIM + m] = U[i][d*DIM + m]*invcount
423 - xav[i*DIM + d]*xav[i*DIM + m];
424 Uaver[3*d+m] += top.atoms.atom[index[i]].m*U[i][d*DIM + m];
427 totmass += top.atoms.atom[index[i]].m;
429 for (d = 0; d < DIM*DIM; d++)
431 Uaver[d] /= totmass;
434 if (bRes)
436 for (d = 0; d < DIM*DIM; d++)
438 average_residues(nullptr, U, d, isize, index, w_rls, &top.atoms);
442 if (bAniso)
444 for (i = 0; i < isize; i++)
446 aid = index[i];
447 pdbatoms->pdbinfo[aid].bAnisotropic = TRUE;
448 pdbatoms->pdbinfo[aid].uij[U11] = static_cast<int>(1e6*U[i][XX*DIM + XX]);
449 pdbatoms->pdbinfo[aid].uij[U22] = static_cast<int>(1e6*U[i][YY*DIM + YY]);
450 pdbatoms->pdbinfo[aid].uij[U33] = static_cast<int>(1e6*U[i][ZZ*DIM + ZZ]);
451 pdbatoms->pdbinfo[aid].uij[U12] = static_cast<int>(1e6*U[i][XX*DIM + YY]);
452 pdbatoms->pdbinfo[aid].uij[U13] = static_cast<int>(1e6*U[i][XX*DIM + ZZ]);
453 pdbatoms->pdbinfo[aid].uij[U23] = static_cast<int>(1e6*U[i][YY*DIM + ZZ]);
456 if (bRes)
458 label = "Residue";
460 else
462 label = "Atom";
465 for (i = 0; i < isize; i++)
467 rmsf[i] = U[i][XX*DIM + XX] + U[i][YY*DIM + YY] + U[i][ZZ*DIM + ZZ];
470 if (dirfn)
472 fprintf(stdout, "\n");
473 print_dir(stdout, Uaver);
474 fp = gmx_ffopen(dirfn, "w");
475 print_dir(fp, Uaver);
476 gmx_ffclose(fp);
479 for (i = 0; i < isize; i++)
481 sfree(U[i]);
483 sfree(U);
485 /* Write RMSF output */
486 if (bReadPDB)
488 bfac = 8.0*M_PI*M_PI/3.0*100;
489 fp = xvgropen(ftp2fn(efXVG, NFILE, fnm), "B-Factors",
490 label, "(A\\b\\S\\So\\N\\S2\\N)", oenv);
491 xvgr_legend(fp, 2, leg, oenv);
492 for (i = 0; (i < isize); i++)
494 if (!bRes || i+1 == isize ||
495 top.atoms.atom[index[i]].resind != top.atoms.atom[index[i+1]].resind)
497 resind = top.atoms.atom[index[i]].resind;
498 pdb_bfac = find_pdb_bfac(pdbatoms, &top.atoms.resinfo[resind],
499 *(top.atoms.atomname[index[i]]));
501 fprintf(fp, "%5d %10.5f %10.5f\n",
502 bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : index[i]+1, rmsf[i]*bfac,
503 pdb_bfac);
506 xvgrclose(fp);
508 else
510 fp = xvgropen(ftp2fn(efXVG, NFILE, fnm), "RMS fluctuation", label, "(nm)", oenv);
511 for (i = 0; i < isize; i++)
513 if (!bRes || i+1 == isize ||
514 top.atoms.atom[index[i]].resind != top.atoms.atom[index[i+1]].resind)
516 fprintf(fp, "%5d %8.4f\n",
517 bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : index[i]+1, std::sqrt(rmsf[i]));
520 xvgrclose(fp);
523 for (i = 0; i < isize; i++)
525 pdbatoms->pdbinfo[index[i]].bfac = 800*M_PI*M_PI/3.0*rmsf[i];
528 if (devfn)
530 for (i = 0; i < isize; i++)
532 rmsf[i] = (rmsd_x[i][XX]+rmsd_x[i][YY]+rmsd_x[i][ZZ])/count;
534 if (bRes)
536 average_residues(rmsf, nullptr, 0, isize, index, w_rls, &top.atoms);
538 /* Write RMSD output */
539 fp = xvgropen(devfn, "RMS Deviation", label, "(nm)", oenv);
540 for (i = 0; i < isize; i++)
542 if (!bRes || i+1 == isize ||
543 top.atoms.atom[index[i]].resind != top.atoms.atom[index[i+1]].resind)
545 fprintf(fp, "%5d %8.4f\n",
546 bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : index[i]+1, std::sqrt(rmsf[i]));
549 xvgrclose(fp);
552 if (opt2bSet("-oq", NFILE, fnm))
554 /* Write a .pdb file with B-factors and optionally anisou records */
555 for (i = 0; i < isize; i++)
557 rvec_inc(xref[index[i]], xcm);
559 write_sto_conf_indexed(opt2fn("-oq", NFILE, fnm), title, pdbatoms, pdbx,
560 nullptr, ePBC, pdbbox, isize, index);
562 if (opt2bSet("-ox", NFILE, fnm))
564 /* Misuse xref as a temporary array */
565 for (i = 0; i < isize; i++)
567 for (d = 0; d < DIM; d++)
569 xref[index[i]][d] = xcm[d] + xav[i*DIM + d];
572 /* Write a .pdb file with B-factors and optionally anisou records */
573 write_sto_conf_indexed(opt2fn("-ox", NFILE, fnm), title, pdbatoms, xref, nullptr,
574 ePBC, pdbbox, isize, index);
576 if (bAniso)
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");
582 if (devfn)
584 do_view(oenv, opt2fn("-od", NFILE, fnm), "-nxy");
587 return 0;