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