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