223a1d7d5995684f3ace3de3cb19a4dfedb0b48d
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_rmsf.c
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, 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 <math.h>
40
41 #include "gromacs/commandline/pargs.h"
42 #include "gromacs/fileio/confio.h"
43 #include "gromacs/fileio/pdbio.h"
44 #include "gromacs/fileio/trxio.h"
45 #include "gromacs/fileio/xvgr.h"
46 #include "gromacs/gmxana/gmx_ana.h"
47 #include "gromacs/gmxana/princ.h"
48 #include "gromacs/legacyheaders/macros.h"
49 #include "gromacs/legacyheaders/typedefs.h"
50 #include "gromacs/legacyheaders/viewit.h"
51 #include "gromacs/linearalgebra/eigensolver.h"
52 #include "gromacs/math/do_fit.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/pbcutil/rmpbc.h"
55 #include "gromacs/topology/index.h"
56 #include "gromacs/utility/cstringutil.h"
57 #include "gromacs/utility/futil.h"
58 #include "gromacs/utility/smalloc.h"
59
60 static real find_pdb_bfac(t_atoms *atoms, t_resinfo *ri, char *atomnm)
61 {
62     char rresnm[8];
63     int  i;
64
65     strcpy(rresnm, *ri->name);
66     rresnm[3] = '\0';
67     for (i = 0; (i < atoms->nr); i++)
68     {
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))
73         {
74             break;
75         }
76     }
77     if (i == atoms->nr)
78     {
79         fprintf(stderr, "\rCan not find %s%d-%s in pdbfile\n",
80                 rresnm, ri->nr, atomnm);
81         return 0.0;
82     }
83
84     return atoms->pdbinfo[i].bfac;
85 }
86
87 void correlate_aniso(const char *fn, t_atoms *ref, t_atoms *calc,
88                      const output_env_t oenv)
89 {
90     FILE *fp;
91     int   i, j;
92
93     fp = xvgropen(fn, "Correlation between X-Ray and Computed Uij", "X-Ray",
94                   "Computed", oenv);
95     for (i = 0; (i < ref->nr); i++)
96     {
97         if (ref->pdbinfo[i].bAnisotropic)
98         {
99             for (j = U11; (j <= U23); j++)
100             {
101                 fprintf(fp, "%10d  %10d\n", ref->pdbinfo[i].uij[j], calc->pdbinfo[i].uij[j]);
102             }
103         }
104     }
105     xvgrclose(fp);
106 }
107
108 static void average_residues(double f[], double **U, int uind,
109                              int isize, atom_id index[], real w_rls[],
110                              t_atoms *atoms)
111 {
112     int    i, j, start;
113     double av, m;
114
115     start = 0;
116     av    = 0;
117     m     = 0;
118     for (i = 0; i < isize; i++)
119     {
120         av += w_rls[index[i]]*(f != NULL ? f[i] : U[i][uind]);
121         m  += w_rls[index[i]];
122         if (i+1 == isize ||
123             atoms->atom[index[i]].resind != atoms->atom[index[i+1]].resind)
124         {
125             av /= m;
126             if (f != NULL)
127             {
128                 for (j = start; j <= i; j++)
129                 {
130                     f[i] = av;
131                 }
132             }
133             else
134             {
135                 for (j = start; j <= i; j++)
136                 {
137                     U[j][uind] = av;
138                 }
139             }
140             start = i+1;
141             av    = 0;
142             m     = 0;
143         }
144     }
145 }
146
147 void print_dir(FILE *fp, real *Uaver)
148 {
149     real eigvec[DIM*DIM];
150     real tmp[DIM*DIM];
151     rvec eigval;
152     int  d, m;
153
154     fprintf(fp, "MSF     X         Y         Z\n");
155     for (d = 0; d < DIM; d++)
156     {
157         fprintf(fp, " %c ", 'X'+d-XX);
158         for (m = 0; m < DIM; m++)
159         {
160             fprintf(fp, " %9.2e", Uaver[3*m+d]);
161         }
162         fprintf(fp, "%s\n", m == DIM ? " (nm^2)" : "");
163     }
164
165     for (m = 0; m < DIM*DIM; m++)
166     {
167         tmp[m] = Uaver[m];
168     }
169
170
171     eigensolver(tmp, DIM, 0, DIM, eigval, eigvec);
172
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++)
177     {
178         fprintf(fp, "  %c   ", 'X'+d-XX);
179         for (m = DIM-1; m >= 0; m--)
180         {
181             fprintf(fp, "%7.4f  ", eigvec[3*m+d]);
182         }
183         fprintf(fp, "\n");
184     }
185 }
186
187 int gmx_rmsf(int argc, char *argv[])
188 {
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 [REF].pdb[ref] file with the coordinates, of the",
195         "structure file, or of a [REF].pdb[ref] file when [TT]-q[tt] is specified.",
196         "Option [TT]-ox[tt] writes the B-factors to a file with the average",
197         "coordinates.[PAR]",
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 [REF].pdb[ref] 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 [REF].pdb[ref] input file is passed to the program and the [TT]-aniso[tt]",
207         "flag is set",
208         "a correlation plot of the Uij will be created, if any anisotropic",
209         "temperature factors are present in the [REF].pdb[ref] 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",
212         "the least."
213     };
214     static gmx_bool  bRes    = FALSE, bAniso = FALSE, bdevX = FALSE, bFit = TRUE;
215     t_pargs          pargs[] = {
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." }
222     };
223     int              natom;
224     int              step, nre, natoms, i, g, m, teller = 0;
225     real             t, lambda, *w_rls, *w_rms;
226
227     t_inputrec       ir;
228     t_topology       top;
229     int              ePBC;
230     t_atoms         *pdbatoms, *refatoms;
231     gmx_bool         bCont;
232
233     matrix           box, pdbbox;
234     rvec            *x, *pdbx, *xref;
235     t_trxstatus     *status;
236     int              npdbatoms, res0;
237     char             buf[256];
238     const char      *label;
239     char             title[STRLEN];
240
241     FILE            *fp;          /* the graphics file */
242     const char      *devfn, *dirfn;
243     int              resind;
244
245     gmx_bool         bReadPDB;
246     atom_id         *index;
247     int              isize;
248     char            *grpnames;
249
250     real             bfac, pdb_bfac, *Uaver;
251     double         **U, *xav;
252     atom_id          aid;
253     rvec            *rmsd_x = NULL;
254     double          *rmsf, invcount, totmass;
255     int              d;
256     real             count = 0;
257     rvec             xcm;
258     gmx_rmpbc_t      gpbc = NULL;
259
260     output_env_t     oenv;
261
262     const char      *leg[2] = { "MD", "X-Ray" };
263
264     t_filenm         fnm[] = {
265         { efTRX, "-f",  NULL,     ffREAD  },
266         { efTPS, NULL,  NULL,     ffREAD  },
267         { efNDX, NULL,  NULL,     ffOPTRD },
268         { efPDB, "-q",  NULL,     ffOPTRD },
269         { efPDB, "-oq", "bfac",   ffOPTWR },
270         { efPDB, "-ox", "xaver",  ffOPTWR },
271         { efXVG, "-o",  "rmsf",   ffWRITE },
272         { efXVG, "-od", "rmsdev", ffOPTWR },
273         { efXVG, "-oc", "correl", ffOPTWR },
274         { efLOG, "-dir", "rmsf",  ffOPTWR }
275     };
276 #define NFILE asize(fnm)
277
278     if (!parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW,
279                            NFILE, fnm, asize(pargs), pargs, asize(desc), desc, 0, NULL,
280                            &oenv))
281     {
282         return 0;
283     }
284
285     bReadPDB = ftp2bSet(efPDB, NFILE, fnm);
286     devfn    = opt2fn_null("-od", NFILE, fnm);
287     dirfn    = opt2fn_null("-dir", NFILE, fnm);
288
289     read_tps_conf(ftp2fn(efTPS, NFILE, fnm), title, &top, &ePBC, &xref, NULL, box, TRUE);
290     snew(w_rls, top.atoms.nr);
291
292     fprintf(stderr, "Select group(s) for root mean square calculation\n");
293     get_index(&top.atoms, ftp2fn_null(efNDX, NFILE, fnm), 1, &isize, &index, &grpnames);
294
295     /* Set the weight */
296     for (i = 0; i < isize; i++)
297     {
298         w_rls[index[i]] = top.atoms.atom[index[i]].m;
299     }
300
301     /* Malloc the rmsf arrays */
302     snew(xav, isize*DIM);
303     snew(U, isize);
304     for (i = 0; i < isize; i++)
305     {
306         snew(U[i], DIM*DIM);
307     }
308     snew(rmsf, isize);
309     if (devfn)
310     {
311         snew(rmsd_x, isize);
312     }
313
314     if (bReadPDB)
315     {
316         get_stx_coordnum(opt2fn("-q", NFILE, fnm), &npdbatoms);
317         snew(pdbatoms, 1);
318         snew(refatoms, 1);
319         init_t_atoms(pdbatoms, npdbatoms, TRUE);
320         init_t_atoms(refatoms, npdbatoms, TRUE);
321         snew(pdbx, npdbatoms);
322         /* Read coordinates twice */
323         read_stx_conf(opt2fn("-q", NFILE, fnm), title, pdbatoms, pdbx, NULL, NULL, pdbbox);
324         read_stx_conf(opt2fn("-q", NFILE, fnm), title, refatoms, pdbx, NULL, NULL, pdbbox);
325     }
326     else
327     {
328         pdbatoms  = &top.atoms;
329         refatoms  = &top.atoms;
330         pdbx      = xref;
331         npdbatoms = pdbatoms->nr;
332         snew(pdbatoms->pdbinfo, npdbatoms);
333         copy_mat(box, pdbbox);
334     }
335
336     if (bFit)
337     {
338         sub_xcm(xref, isize, index, top.atoms.atom, xcm, FALSE);
339     }
340
341     natom = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
342
343     if (bFit)
344     {
345         gpbc = gmx_rmpbc_init(&top.idef, ePBC, natom);
346     }
347
348     /* Now read the trj again to compute fluctuations */
349     teller = 0;
350     do
351     {
352         if (bFit)
353         {
354             /* Remove periodic boundary */
355             gmx_rmpbc(gpbc, natom, box, x);
356
357             /* Set center of mass to zero */
358             sub_xcm(x, isize, index, top.atoms.atom, xcm, FALSE);
359
360             /* Fit to reference structure */
361             do_fit(natom, w_rls, xref, x);
362         }
363
364         /* Calculate Anisotropic U Tensor */
365         for (i = 0; i < isize; i++)
366         {
367             aid = index[i];
368             for (d = 0; d < DIM; d++)
369             {
370                 xav[i*DIM + d] += x[aid][d];
371                 for (m = 0; m < DIM; m++)
372                 {
373                     U[i][d*DIM + m] += x[aid][d]*x[aid][m];
374                 }
375             }
376         }
377
378         if (devfn)
379         {
380             /* Calculate RMS Deviation */
381             for (i = 0; (i < isize); i++)
382             {
383                 aid = index[i];
384                 for (d = 0; (d < DIM); d++)
385                 {
386                     rmsd_x[i][d] += sqr(x[aid][d]-xref[aid][d]);
387                 }
388             }
389         }
390         count += 1.0;
391         teller++;
392     }
393     while (read_next_x(oenv, status, &t, x, box));
394     close_trj(status);
395
396     if (bFit)
397     {
398         gmx_rmpbc_done(gpbc);
399     }
400
401
402     invcount = 1.0/count;
403     snew(Uaver, DIM*DIM);
404     totmass = 0;
405     for (i = 0; i < isize; i++)
406     {
407         for (d = 0; d < DIM; d++)
408         {
409             xav[i*DIM + d] *= invcount;
410         }
411         for (d = 0; d < DIM; d++)
412         {
413             for (m = 0; m < DIM; m++)
414             {
415                 U[i][d*DIM + m] = U[i][d*DIM + m]*invcount
416                     - xav[i*DIM + d]*xav[i*DIM + m];
417                 Uaver[3*d+m] += top.atoms.atom[index[i]].m*U[i][d*DIM + m];
418             }
419         }
420         totmass += top.atoms.atom[index[i]].m;
421     }
422     for (d = 0; d < DIM*DIM; d++)
423     {
424         Uaver[d] /= totmass;
425     }
426
427     if (bRes)
428     {
429         for (d = 0; d < DIM*DIM; d++)
430         {
431             average_residues(NULL, U, d, isize, index, w_rls, &top.atoms);
432         }
433     }
434
435     if (bAniso)
436     {
437         for (i = 0; i < isize; i++)
438         {
439             aid = index[i];
440             pdbatoms->pdbinfo[aid].bAnisotropic = TRUE;
441             pdbatoms->pdbinfo[aid].uij[U11]     = 1e6*U[i][XX*DIM + XX];
442             pdbatoms->pdbinfo[aid].uij[U22]     = 1e6*U[i][YY*DIM + YY];
443             pdbatoms->pdbinfo[aid].uij[U33]     = 1e6*U[i][ZZ*DIM + ZZ];
444             pdbatoms->pdbinfo[aid].uij[U12]     = 1e6*U[i][XX*DIM + YY];
445             pdbatoms->pdbinfo[aid].uij[U13]     = 1e6*U[i][XX*DIM + ZZ];
446             pdbatoms->pdbinfo[aid].uij[U23]     = 1e6*U[i][YY*DIM + ZZ];
447         }
448     }
449     if (bRes)
450     {
451         label = "Residue";
452     }
453     else
454     {
455         label = "Atom";
456     }
457
458     for (i = 0; i < isize; i++)
459     {
460         rmsf[i] = U[i][XX*DIM + XX] + U[i][YY*DIM + YY] + U[i][ZZ*DIM + ZZ];
461     }
462
463     if (dirfn)
464     {
465         fprintf(stdout, "\n");
466         print_dir(stdout, Uaver);
467         fp = gmx_ffopen(dirfn, "w");
468         print_dir(fp, Uaver);
469         gmx_ffclose(fp);
470     }
471
472     for (i = 0; i < isize; i++)
473     {
474         sfree(U[i]);
475     }
476     sfree(U);
477
478     /* Write RMSF output */
479     if (bReadPDB)
480     {
481         bfac = 8.0*M_PI*M_PI/3.0*100;
482         fp   = xvgropen(ftp2fn(efXVG, NFILE, fnm), "B-Factors",
483                         label, "(A\\b\\S\\So\\N\\S2\\N)", oenv);
484         xvgr_legend(fp, 2, leg, oenv);
485         for (i = 0; (i < isize); i++)
486         {
487             if (!bRes || i+1 == isize ||
488                 top.atoms.atom[index[i]].resind != top.atoms.atom[index[i+1]].resind)
489             {
490                 resind    = top.atoms.atom[index[i]].resind;
491                 pdb_bfac  = find_pdb_bfac(pdbatoms, &top.atoms.resinfo[resind],
492                                           *(top.atoms.atomname[index[i]]));
493
494                 fprintf(fp, "%5d  %10.5f  %10.5f\n",
495                         bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : index[i]+1, rmsf[i]*bfac,
496                         pdb_bfac);
497             }
498         }
499         xvgrclose(fp);
500     }
501     else
502     {
503         fp = xvgropen(ftp2fn(efXVG, NFILE, fnm), "RMS fluctuation", label, "(nm)", oenv);
504         for (i = 0; i < isize; i++)
505         {
506             if (!bRes || i+1 == isize ||
507                 top.atoms.atom[index[i]].resind != top.atoms.atom[index[i+1]].resind)
508             {
509                 fprintf(fp, "%5d %8.4f\n",
510                         bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : index[i]+1, sqrt(rmsf[i]));
511             }
512         }
513         xvgrclose(fp);
514     }
515
516     for (i = 0; i < isize; i++)
517     {
518         pdbatoms->pdbinfo[index[i]].bfac = 800*M_PI*M_PI/3.0*rmsf[i];
519     }
520
521     if (devfn)
522     {
523         for (i = 0; i < isize; i++)
524         {
525             rmsf[i] = (rmsd_x[i][XX]+rmsd_x[i][YY]+rmsd_x[i][ZZ])/count;
526         }
527         if (bRes)
528         {
529             average_residues(rmsf, NULL, 0, isize, index, w_rls, &top.atoms);
530         }
531         /* Write RMSD output */
532         fp = xvgropen(devfn, "RMS Deviation", label, "(nm)", oenv);
533         for (i = 0; i < isize; i++)
534         {
535             if (!bRes || i+1 == isize ||
536                 top.atoms.atom[index[i]].resind != top.atoms.atom[index[i+1]].resind)
537             {
538                 fprintf(fp, "%5d %8.4f\n",
539                         bRes ? top.atoms.resinfo[top.atoms.atom[index[i]].resind].nr : index[i]+1, sqrt(rmsf[i]));
540             }
541         }
542         xvgrclose(fp);
543     }
544
545     if (opt2bSet("-oq", NFILE, fnm))
546     {
547         /* Write a .pdb file with B-factors and optionally anisou records */
548         for (i = 0; i < isize; i++)
549         {
550             rvec_inc(xref[index[i]], xcm);
551         }
552         write_sto_conf_indexed(opt2fn("-oq", NFILE, fnm), title, pdbatoms, pdbx,
553                                NULL, ePBC, pdbbox, isize, index);
554     }
555     if (opt2bSet("-ox", NFILE, fnm))
556     {
557         /* Misuse xref as a temporary array */
558         for (i = 0; i < isize; i++)
559         {
560             for (d = 0; d < DIM; d++)
561             {
562                 xref[index[i]][d] = xcm[d] + xav[i*DIM + d];
563             }
564         }
565         /* Write a .pdb file with B-factors and optionally anisou records */
566         write_sto_conf_indexed(opt2fn("-ox", NFILE, fnm), title, pdbatoms, xref, NULL,
567                                ePBC, pdbbox, isize, index);
568     }
569     if (bAniso)
570     {
571         correlate_aniso(opt2fn("-oc", NFILE, fnm), refatoms, pdbatoms, oenv);
572         do_view(oenv, opt2fn("-oc", NFILE, fnm), "-nxy");
573     }
574     do_view(oenv, opt2fn("-o", NFILE, fnm), "-nxy");
575     if (devfn)
576     {
577         do_view(oenv, opt2fn("-od", NFILE, fnm), "-nxy");
578     }
579
580     return 0;
581 }