Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxana / gmx_rmsf.c
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  *
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  *
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  *
30  * For more info, check our website at http://www.gromacs.org
31  *
32  * And Hey:
33  * Green Red Orange Magenta Azure Cyan Skyblue
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include "smalloc.h"
40 #include <math.h>
41 #include "macros.h"
42 #include "typedefs.h"
43 #include "xvgr.h"
44 #include "copyrite.h"
45 #include "statutil.h"
46 #include "string2.h"
47 #include "vec.h"
48 #include "index.h"
49 #include "pdbio.h"
50 #include "tpxio.h"
51 #include "pbc.h"
52 #include "gmx_fatal.h"
53 #include "futil.h"
54 #include "do_fit.h"
55 #include "princ.h"
56 #include "rmpbc.h"
57 #include "confio.h"
58 #include "gmx_ana.h"
59
60 #include "gromacs/linearalgebra/eigensolver.h"
61
62 static real find_pdb_bfac(t_atoms *atoms, t_resinfo *ri, char *atomnm)
63 {
64     char rresnm[8];
65     int  i;
66
67     strcpy(rresnm, *ri->name);
68     rresnm[3] = '\0';
69     for (i = 0; (i < atoms->nr); i++)
70     {
71         if ((ri->nr == atoms->resinfo[atoms->atom[i].resind].nr) &&
72             (ri->ic == atoms->resinfo[atoms->atom[i].resind].ic) &&
73             (strcmp(*atoms->resinfo[atoms->atom[i].resind].name, rresnm) == 0) &&
74             (strstr(*atoms->atomname[i], atomnm) != NULL))
75         {
76             break;
77         }
78     }
79     if (i == atoms->nr)
80     {
81         fprintf(stderr, "\rCan not find %s%d-%s in pdbfile\n",
82                 rresnm, ri->nr, atomnm);
83         return 0.0;
84     }
85
86     return atoms->pdbinfo[i].bfac;
87 }
88
89 void correlate_aniso(const char *fn, t_atoms *ref, t_atoms *calc,
90                      const output_env_t oenv)
91 {
92     FILE *fp;
93     int   i, j;
94
95     fp = xvgropen(fn, "Correlation between X-Ray and Computed Uij", "X-Ray",
96                   "Computed", oenv);
97     for (i = 0; (i < ref->nr); i++)
98     {
99         if (ref->pdbinfo[i].bAnisotropic)
100         {
101             for (j = U11; (j <= U23); j++)
102             {
103                 fprintf(fp, "%10d  %10d\n", ref->pdbinfo[i].uij[j], calc->pdbinfo[i].uij[j]);
104             }
105         }
106     }
107     ffclose(fp);
108 }
109
110 static void average_residues(double f[], double **U, int uind,
111                              int isize, atom_id index[], real w_rls[],
112                              t_atoms *atoms)
113 {
114     int    i, j, start;
115     double av, m;
116
117     start = 0;
118     av    = 0;
119     m     = 0;
120     for (i = 0; i < isize; i++)
121     {
122         av += w_rls[index[i]]*(f != NULL ? f[i] : U[i][uind]);
123         m  += w_rls[index[i]];
124         if (i+1 == isize ||
125             atoms->atom[index[i]].resind != atoms->atom[index[i+1]].resind)
126         {
127             av /= m;
128             if (f != NULL)
129             {
130                 for (j = start; j <= i; j++)
131                 {
132                     f[i] = av;
133                 }
134             }
135             else
136             {
137                 for (j = start; j <= i; j++)
138                 {
139                     U[j][uind] = av;
140                 }
141             }
142             start = i+1;
143             av    = 0;
144             m     = 0;
145         }
146     }
147 }
148
149 void print_dir(FILE *fp, real *Uaver)
150 {
151     real eigvec[DIM*DIM];
152     real tmp[DIM*DIM];
153     rvec eigval;
154     int  d, m;
155
156     fprintf(fp, "MSF     X         Y         Z\n");
157     for (d = 0; d < DIM; d++)
158     {
159         fprintf(fp, " %c ", 'X'+d-XX);
160         for (m = 0; m < DIM; m++)
161         {
162             fprintf(fp, " %9.2e", Uaver[3*m+d]);
163         }
164         fprintf(fp, "%s\n", m == DIM ? " (nm^2)" : "");
165     }
166
167     for (m = 0; m < DIM*DIM; m++)
168     {
169         tmp[m] = Uaver[m];
170     }
171
172
173     eigensolver(tmp, DIM, 0, DIM, eigval, eigvec);
174
175     fprintf(fp, "\n             Eigenvectors\n\n");
176     fprintf(fp, "Eigv  %-8.2e %-8.2e %-8.2e (nm^2)\n\n",
177             eigval[2], eigval[1], eigval[0]);
178     for (d = 0; d < DIM; d++)
179     {
180         fprintf(fp, "  %c   ", 'X'+d-XX);
181         for (m = DIM-1; m >= 0; m--)
182         {
183             fprintf(fp, "%7.4f  ", eigvec[3*m+d]);
184         }
185         fprintf(fp, "\n");
186     }
187 }
188
189 int gmx_rmsf(int argc, char *argv[])
190 {
191     const char      *desc[] = {
192         "[TT]g_rmsf[tt] computes the root mean square fluctuation (RMSF, i.e. standard ",
193         "deviation) of atomic positions in the trajectory (supplied with [TT]-f[tt])",
194         "after (optionally) fitting to a reference frame (supplied with [TT]-s[tt]).[PAR]",
195         "With option [TT]-oq[tt] the RMSF values are converted to B-factor",
196         "values, which are written to a [TT].pdb[tt] file with the coordinates, of the",
197         "structure file, or of a [TT].pdb[tt] file when [TT]-q[tt] is specified.",
198         "Option [TT]-ox[tt] writes the B-factors to a file with the average",
199         "coordinates.[PAR]",
200         "With the option [TT]-od[tt] the root mean square deviation with",
201         "respect to the reference structure is calculated.[PAR]",
202         "With the option [TT]-aniso[tt], [TT]g_rmsf[tt] will compute anisotropic",
203         "temperature factors and then it will also output average coordinates",
204         "and a [TT].pdb[tt] file with ANISOU records (corresonding to the [TT]-oq[tt]",
205         "or [TT]-ox[tt] option). Please note that the U values",
206         "are orientation-dependent, so before comparison with experimental data",
207         "you should verify that you fit to the experimental coordinates.[PAR]",
208         "When a [TT].pdb[tt] input file is passed to the program and the [TT]-aniso[tt]",
209         "flag is set",
210         "a correlation plot of the Uij will be created, if any anisotropic",
211         "temperature factors are present in the [TT].pdb[tt] file.[PAR]",
212         "With option [TT]-dir[tt] the average MSF (3x3) matrix is diagonalized.",
213         "This shows the directions in which the atoms fluctuate the most and",
214         "the least."
215     };
216     static gmx_bool  bRes    = FALSE, bAniso = FALSE, bdevX = FALSE, bFit = TRUE;
217     t_pargs          pargs[] = {
218         { "-res", FALSE, etBOOL, {&bRes},
219           "Calculate averages for each residue" },
220         { "-aniso", FALSE, etBOOL, {&bAniso},
221           "Compute anisotropic termperature factors" },
222         { "-fit", FALSE, etBOOL, {&bFit},
223           "Do a least squares superposition before computing RMSF. Without this you must make sure that the reference structure and the trajectory match." }
224     };
225     int              natom;
226     int              step, nre, natoms, i, g, m, teller = 0;
227     real             t, lambda, *w_rls, *w_rms;
228
229     t_tpxheader      header;
230     t_inputrec       ir;
231     t_topology       top;
232     int              ePBC;
233     t_atoms         *pdbatoms, *refatoms;
234     gmx_bool         bCont;
235
236     matrix           box, pdbbox;
237     rvec            *x, *pdbx, *xref;
238     t_trxstatus     *status;
239     int              npdbatoms, res0;
240     char             buf[256];
241     const char      *label;
242     char             title[STRLEN];
243
244     FILE            *fp;          /* the graphics file */
245     const char      *devfn, *dirfn;
246     int              resind;
247
248     gmx_bool         bReadPDB;
249     atom_id         *index;
250     int              isize;
251     char            *grpnames;
252
253     real             bfac, pdb_bfac, *Uaver;
254     double         **U, *xav;
255     atom_id          aid;
256     rvec            *rmsd_x = NULL;
257     double          *rmsf, invcount, totmass;
258     int              d;
259     real             count = 0;
260     rvec             xcm;
261     gmx_rmpbc_t      gpbc = NULL;
262
263     output_env_t     oenv;
264
265     const char      *leg[2] = { "MD", "X-Ray" };
266
267     t_filenm         fnm[] = {
268         { efTRX, "-f",  NULL,     ffREAD  },
269         { efTPS, NULL,  NULL,     ffREAD  },
270         { efNDX, NULL,  NULL,     ffOPTRD },
271         { efPDB, "-q",  NULL,     ffOPTRD },
272         { efPDB, "-oq", "bfac",   ffOPTWR },
273         { efPDB, "-ox", "xaver",  ffOPTWR },
274         { efXVG, "-o",  "rmsf",   ffWRITE },
275         { efXVG, "-od", "rmsdev", ffOPTWR },
276         { efXVG, "-oc", "correl", ffOPTWR },
277         { efLOG, "-dir", "rmsf",  ffOPTWR }
278     };
279 #define NFILE asize(fnm)
280
281     parse_common_args(&argc, argv, PCA_CAN_TIME | PCA_CAN_VIEW | PCA_BE_NICE,
282                       NFILE, fnm, asize(pargs), pargs, asize(desc), desc, 0, NULL,
283                       &oenv);
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, box);
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, natom, 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 = ffopen(dirfn, "w");
468         print_dir(fp, Uaver);
469         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         ffclose(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         ffclose(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         ffclose(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     thanx(stderr);
581
582     return 0;
583 }