354f092b02a9539586ee7cd811ede5568450044f
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_rmsf.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include <cassert>
40 #include <cmath>
41 #include <cstring>
42
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"
63
64 static real find_pdb_bfac(const t_atoms *atoms, t_resinfo *ri, char *atomnm)
65 {
66     char rresnm[8];
67     int  i;
68
69     std::strcpy(rresnm, *ri->name);
70     rresnm[3] = '\0';
71     for (i = 0; (i < atoms->nr); i++)
72     {
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))
77         {
78             break;
79         }
80     }
81     if (i == atoms->nr)
82     {
83         fprintf(stderr, "\rCan not find %s%d-%s in pdbfile\n",
84                 rresnm, ri->nr, atomnm);
85         fflush(stderr);
86         return 0.0;
87     }
88
89     return atoms->pdbinfo[i].bfac;
90 }
91
92 void correlate_aniso(const char *fn, t_atoms *ref, t_atoms *calc,
93                      const gmx_output_env_t *oenv)
94 {
95     FILE *fp;
96     int   i, j;
97
98     fp = xvgropen(fn, "Correlation between X-Ray and Computed Uij", "X-Ray",
99                   "Computed", oenv);
100     for (i = 0; (i < ref->nr); i++)
101     {
102         if (ref->pdbinfo[i].bAnisotropic)
103         {
104             for (j = U11; (j <= U23); j++)
105             {
106                 fprintf(fp, "%10d  %10d\n", ref->pdbinfo[i].uij[j], calc->pdbinfo[i].uij[j]);
107             }
108         }
109     }
110     xvgrclose(fp);
111 }
112
113 static void average_residues(double f[], double **U, int uind,
114                              int isize, int index[], real w_rls[],
115                              const t_atoms *atoms)
116 {
117     int    i, j, start;
118     double av, m;
119
120     start = 0;
121     av    = 0;
122     m     = 0;
123     if (!f)
124     {
125         assert(U);
126     }
127     for (i = 0; i < isize; i++)
128     {
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)
133         {
134             av /= m;
135             if (f != nullptr)
136             {
137                 for (j = start; j <= i; j++)
138                 {
139                     f[i] = av;
140                 }
141             }
142             else
143             {
144                 for (j = start; j <= i; j++)
145                 {
146                     U[j][uind] = av;
147                 }
148             }
149             start = i+1;
150             av    = 0;
151             m     = 0;
152         }
153     }
154 }
155
156 void print_dir(FILE *fp, real *Uaver)
157 {
158     real eigvec[DIM*DIM];
159     real tmp[DIM*DIM];
160     rvec eigval;
161     int  d, m;
162
163     fprintf(fp, "MSF     X         Y         Z\n");
164     for (d = 0; d < DIM; d++)
165     {
166         fprintf(fp, " %c ", 'X'+d-XX);
167         for (m = 0; m < DIM; m++)
168         {
169             fprintf(fp, " %9.2e", Uaver[3*m+d]);
170         }
171         fprintf(fp, "%s\n", m == DIM ? " (nm^2)" : "");
172     }
173
174     for (m = 0; m < DIM*DIM; m++)
175     {
176         tmp[m] = Uaver[m];
177     }
178
179
180     eigensolver(tmp, DIM, 0, DIM, eigval, eigvec);
181
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++)
186     {
187         fprintf(fp, "  %c   ", 'X'+d-XX);
188         for (m = DIM-1; m >= 0; m--)
189         {
190             fprintf(fp, "%7.4f  ", eigvec[3*m+d]);
191         }
192         fprintf(fp, "\n");
193     }
194 }
195
196 int gmx_rmsf(int argc, char *argv[])
197 {
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."
222     };
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." }
231     };
232     int               natom;
233     int               i, m, teller = 0;
234     real              t, *w_rls;
235
236     t_topology        top;
237     int               ePBC;
238     t_atoms          *pdbatoms, *refatoms;
239
240     matrix            box, pdbbox;
241     rvec             *x, *pdbx, *xref;
242     t_trxstatus      *status;
243     const char       *label;
244
245     FILE             *fp;         /* the graphics file */
246     const char       *devfn, *dirfn;
247     int               resind;
248
249     gmx_bool          bReadPDB;
250     int              *index;
251     int               isize;
252     char             *grpnames;
253
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;
263
264     gmx_output_env_t *oenv;
265
266     const char       *leg[2] = { "MD", "X-Ray" };
267
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 }
279     };
280 #define NFILE asize(fnm)
281
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))
285     {
286         return 0;
287     }
288
289     bReadPDB = ftp2bSet(efPDB, NFILE, fnm);
290     devfn    = opt2fn_null("-od", NFILE, fnm);
291     dirfn    = opt2fn_null("-dir", NFILE, fnm);
292
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);
296
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);
299
300     /* Set the weight */
301     for (i = 0; i < isize; i++)
302     {
303         w_rls[index[i]] = top.atoms.atom[index[i]].m;
304     }
305
306     /* Malloc the rmsf arrays */
307     snew(xav, isize*DIM);
308     snew(U, isize);
309     for (i = 0; i < isize; i++)
310     {
311         snew(U[i], DIM*DIM);
312     }
313     snew(rmsf, isize);
314     if (devfn)
315     {
316         snew(rmsd_x, isize);
317     }
318
319     if (bReadPDB)
320     {
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         /* TODO Should this assert that top_pdb->atoms.nr == top.atoms.nr?
329          * See discussion at https://gerrit.gromacs.org/#/c/6430/1 */
330         title = *top_pdb->name;
331         snew(refatoms, 1);
332         *refatoms = top_pdb->atoms;
333         sfree(top_pdb);
334     }
335     else
336     {
337         pdbatoms  = &top.atoms;
338         refatoms  = &top.atoms;
339         pdbx      = xref;
340         snew(pdbatoms->pdbinfo, pdbatoms->nr);
341         pdbatoms->havePdbInfo = TRUE;
342         copy_mat(box, pdbbox);
343     }
344
345     if (bFit)
346     {
347         sub_xcm(xref, isize, index, top.atoms.atom, xcm, FALSE);
348     }
349
350     natom = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
351
352     if (bFit)
353     {
354         gpbc = gmx_rmpbc_init(&top.idef, ePBC, natom);
355     }
356
357     /* Now read the trj again to compute fluctuations */
358     teller = 0;
359     do
360     {
361         if (bFit)
362         {
363             /* Remove periodic boundary */
364             gmx_rmpbc(gpbc, natom, box, x);
365
366             /* Set center of mass to zero */
367             sub_xcm(x, isize, index, top.atoms.atom, xcm, FALSE);
368
369             /* Fit to reference structure */
370             do_fit(natom, w_rls, xref, x);
371         }
372
373         /* Calculate Anisotropic U Tensor */
374         for (i = 0; i < isize; i++)
375         {
376             aid = index[i];
377             for (d = 0; d < DIM; d++)
378             {
379                 xav[i*DIM + d] += x[aid][d];
380                 for (m = 0; m < DIM; m++)
381                 {
382                     U[i][d*DIM + m] += x[aid][d]*x[aid][m];
383                 }
384             }
385         }
386
387         if (devfn)
388         {
389             /* Calculate RMS Deviation */
390             for (i = 0; (i < isize); i++)
391             {
392                 aid = index[i];
393                 for (d = 0; (d < DIM); d++)
394                 {
395                     rmsd_x[i][d] += gmx::square(x[aid][d]-xref[aid][d]);
396                 }
397             }
398         }
399         count += 1.0;
400         teller++;
401     }
402     while (read_next_x(oenv, status, &t, x, box));
403     close_trx(status);
404
405     if (bFit)
406     {
407         gmx_rmpbc_done(gpbc);
408     }
409
410
411     invcount = 1.0/count;
412     snew(Uaver, DIM*DIM);
413     totmass = 0;
414     for (i = 0; i < isize; i++)
415     {
416         for (d = 0; d < DIM; d++)
417         {
418             xav[i*DIM + d] *= invcount;
419         }
420         for (d = 0; d < DIM; d++)
421         {
422             for (m = 0; m < DIM; m++)
423             {
424                 U[i][d*DIM + m] = U[i][d*DIM + m]*invcount
425                     - xav[i*DIM + d]*xav[i*DIM + m];
426                 Uaver[3*d+m] += top.atoms.atom[index[i]].m*U[i][d*DIM + m];
427             }
428         }
429         totmass += top.atoms.atom[index[i]].m;
430     }
431     for (d = 0; d < DIM*DIM; d++)
432     {
433         Uaver[d] /= totmass;
434     }
435
436     if (bRes)
437     {
438         for (d = 0; d < DIM*DIM; d++)
439         {
440             average_residues(nullptr, U, d, isize, index, w_rls, &top.atoms);
441         }
442     }
443
444     if (bAniso)
445     {
446         for (i = 0; i < isize; i++)
447         {
448             aid = index[i];
449             pdbatoms->pdbinfo[aid].bAnisotropic = TRUE;
450             pdbatoms->pdbinfo[aid].uij[U11]     = static_cast<int>(1e6*U[i][XX*DIM + XX]);
451             pdbatoms->pdbinfo[aid].uij[U22]     = static_cast<int>(1e6*U[i][YY*DIM + YY]);
452             pdbatoms->pdbinfo[aid].uij[U33]     = static_cast<int>(1e6*U[i][ZZ*DIM + ZZ]);
453             pdbatoms->pdbinfo[aid].uij[U12]     = static_cast<int>(1e6*U[i][XX*DIM + YY]);
454             pdbatoms->pdbinfo[aid].uij[U13]     = static_cast<int>(1e6*U[i][XX*DIM + ZZ]);
455             pdbatoms->pdbinfo[aid].uij[U23]     = static_cast<int>(1e6*U[i][YY*DIM + ZZ]);
456         }
457     }
458     if (bRes)
459     {
460         label = "Residue";
461     }
462     else
463     {
464         label = "Atom";
465     }
466
467     for (i = 0; i < isize; i++)
468     {
469         rmsf[i] = U[i][XX*DIM + XX] + U[i][YY*DIM + YY] + U[i][ZZ*DIM + ZZ];
470     }
471
472     if (dirfn)
473     {
474         fprintf(stdout, "\n");
475         print_dir(stdout, Uaver);
476         fp = gmx_ffopen(dirfn, "w");
477         print_dir(fp, Uaver);
478         gmx_ffclose(fp);
479     }
480
481     for (i = 0; i < isize; i++)
482     {
483         sfree(U[i]);
484     }
485     sfree(U);
486
487     /* Write RMSF output */
488     if (bReadPDB)
489     {
490         bfac = 8.0*M_PI*M_PI/3.0*100;
491         fp   = xvgropen(ftp2fn(efXVG, NFILE, fnm), "B-Factors",
492                         label, "(A\\b\\S\\So\\N\\S2\\N)", oenv);
493         xvgr_legend(fp, 2, leg, oenv);
494         for (i = 0; (i < isize); i++)
495         {
496             if (!bRes || i+1 == isize ||
497                 top.atoms.atom[index[i]].resind != top.atoms.atom[index[i+1]].resind)
498             {
499                 resind    = top.atoms.atom[index[i]].resind;
500                 pdb_bfac  = find_pdb_bfac(pdbatoms, &top.atoms.resinfo[resind],
501                                           *(top.atoms.atomname[index[i]]));
502
503                 fprintf(fp, "%5d  %10.5f  %10.5f\n",
504                         bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : index[i]+1, rmsf[i]*bfac,
505                         pdb_bfac);
506             }
507         }
508         xvgrclose(fp);
509     }
510     else
511     {
512         fp = xvgropen(ftp2fn(efXVG, NFILE, fnm), "RMS fluctuation", label, "(nm)", oenv);
513         for (i = 0; i < isize; i++)
514         {
515             if (!bRes || i+1 == isize ||
516                 top.atoms.atom[index[i]].resind != top.atoms.atom[index[i+1]].resind)
517             {
518                 fprintf(fp, "%5d %8.4f\n",
519                         bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : index[i]+1, std::sqrt(rmsf[i]));
520             }
521         }
522         xvgrclose(fp);
523     }
524
525     for (i = 0; i < isize; i++)
526     {
527         pdbatoms->pdbinfo[index[i]].bfac = 800*M_PI*M_PI/3.0*rmsf[i];
528     }
529
530     if (devfn)
531     {
532         for (i = 0; i < isize; i++)
533         {
534             rmsf[i] = (rmsd_x[i][XX]+rmsd_x[i][YY]+rmsd_x[i][ZZ])/count;
535         }
536         if (bRes)
537         {
538             average_residues(rmsf, nullptr, 0, isize, index, w_rls, &top.atoms);
539         }
540         /* Write RMSD output */
541         fp = xvgropen(devfn, "RMS Deviation", label, "(nm)", oenv);
542         for (i = 0; i < isize; i++)
543         {
544             if (!bRes || i+1 == isize ||
545                 top.atoms.atom[index[i]].resind != top.atoms.atom[index[i+1]].resind)
546             {
547                 fprintf(fp, "%5d %8.4f\n",
548                         bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : index[i]+1, std::sqrt(rmsf[i]));
549             }
550         }
551         xvgrclose(fp);
552     }
553
554     if (opt2bSet("-oq", NFILE, fnm))
555     {
556         /* Write a .pdb file with B-factors and optionally anisou records */
557         for (i = 0; i < isize; i++)
558         {
559             rvec_inc(pdbx[index[i]], xcm);
560         }
561         write_sto_conf_indexed(opt2fn("-oq", NFILE, fnm), title, pdbatoms, pdbx,
562                                nullptr, ePBC, pdbbox, isize, index);
563     }
564     if (opt2bSet("-ox", NFILE, fnm))
565     {
566         rvec *bFactorX;
567         snew(bFactorX, top.atoms.nr);
568         for (i = 0; i < isize; i++)
569         {
570             for (d = 0; d < DIM; d++)
571             {
572                 bFactorX[index[i]][d] = xcm[d] + xav[i*DIM + d];
573             }
574         }
575         /* Write a .pdb file with B-factors and optionally anisou records */
576         write_sto_conf_indexed(opt2fn("-ox", NFILE, fnm), title, pdbatoms, bFactorX, nullptr,
577                                ePBC, pdbbox, isize, index);
578         sfree(bFactorX);
579     }
580     if (bAniso)
581     {
582         correlate_aniso(opt2fn("-oc", NFILE, fnm), refatoms, pdbatoms, oenv);
583         do_view(oenv, opt2fn("-oc", NFILE, fnm), "-nxy");
584     }
585     do_view(oenv, opt2fn("-o", NFILE, fnm), "-nxy");
586     if (devfn)
587     {
588         do_view(oenv, opt2fn("-od", NFILE, fnm), "-nxy");
589     }
590
591     return 0;
592 }