Merge branch 'release-4-6' into master
[alexxy/gromacs.git] / src / gromacs / gmxana / dlist.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 <stdlib.h>
40
41 #include "string2.h"
42 #include "smalloc.h"
43 #include "gstat.h"
44 #include "gmx_fatal.h"
45 #include "index.h"
46
47 t_dlist *mk_dlist(FILE *log,
48                   t_atoms *atoms, int *nlist,
49                   gmx_bool bPhi, gmx_bool bPsi, gmx_bool bChi, gmx_bool bHChi,
50                   int maxchi, int r0, gmx_residuetype_t rt)
51 {
52     int       ires, i, j, k, ii;
53     t_dihatms atm, prev;
54     int       nl = 0, nc[edMax];
55     char     *thisres;
56     t_dlist  *dl;
57
58     snew(dl, atoms->nres+1);
59     prev.C = prev.Cn[1] = -1; /* Keep the compiler quiet */
60     for (i = 0; (i < edMax); i++)
61     {
62         nc[i] = 0;
63     }
64     ires = -1;
65     i    =  0;
66     while (i < atoms->nr)
67     {
68         ires = atoms->atom[i].resind;
69
70         /* Initiate all atom numbers to -1 */
71         atm.minC = atm.H = atm.N = atm.C = atm.O = atm.minCalpha = -1;
72         for (j = 0; (j < MAXCHI+3); j++)
73         {
74             atm.Cn[j] = -1;
75         }
76
77         /* Look for atoms in this residue */
78         /* maybe should allow for chis to hydrogens? */
79         while ((i < atoms->nr) && (atoms->atom[i].resind == ires))
80         {
81             if ((strcmp(*(atoms->atomname[i]), "H") == 0) ||
82                 (strcmp(*(atoms->atomname[i]), "H1") == 0) ||
83                 (strcmp(*(atoms->atomname[i]), "HN") == 0) )
84             {
85                 atm.H = i;
86             }
87             else if (strcmp(*(atoms->atomname[i]), "N") == 0)
88             {
89                 atm.N = i;
90             }
91             else if (strcmp(*(atoms->atomname[i]), "C") == 0)
92             {
93                 atm.C = i;
94             }
95             else if ((strcmp(*(atoms->atomname[i]), "O") == 0) ||
96                      (strcmp(*(atoms->atomname[i]), "O1") == 0))
97             {
98                 atm.O = i;
99             }
100             else if (strcmp(*(atoms->atomname[i]), "CA") == 0)
101             {
102                 atm.Cn[1] = i;
103             }
104             else if (strcmp(*(atoms->atomname[i]), "CB") == 0)
105             {
106                 atm.Cn[2] = i;
107             }
108             else if ((strcmp(*(atoms->atomname[i]), "CG") == 0)  ||
109                      (strcmp(*(atoms->atomname[i]), "CG1") == 0) ||
110                      (strcmp(*(atoms->atomname[i]), "OG") == 0)  ||
111                      (strcmp(*(atoms->atomname[i]), "OG1") == 0) ||
112                      (strcmp(*(atoms->atomname[i]), "SG") == 0))
113             {
114                 atm.Cn[3] = i;
115             }
116             else if ((strcmp(*(atoms->atomname[i]), "CD") == 0)  ||
117                      (strcmp(*(atoms->atomname[i]), "CD1") == 0) ||
118                      (strcmp(*(atoms->atomname[i]), "SD") == 0)  ||
119                      (strcmp(*(atoms->atomname[i]), "OD1") == 0) ||
120                      (strcmp(*(atoms->atomname[i]), "ND1") == 0))
121             {
122                 atm.Cn[4] = i;
123             }
124             /* by grs - split the Cn[4] into 2 bits to check allowing dih to H */
125             else if (bHChi && ((strcmp(*(atoms->atomname[i]), "HG")  == 0) ||
126                                (strcmp(*(atoms->atomname[i]), "HG1")  == 0)) )
127             {
128                 atm.Cn[4] = i;
129             }
130             else if ((strcmp(*(atoms->atomname[i]), "CE") == 0) ||
131                      (strcmp(*(atoms->atomname[i]), "CE1") == 0) ||
132                      (strcmp(*(atoms->atomname[i]), "OE1") == 0) ||
133                      (strcmp(*(atoms->atomname[i]), "NE") == 0))
134             {
135                 atm.Cn[5] = i;
136             }
137             else if ((strcmp(*(atoms->atomname[i]), "CZ") == 0) ||
138                      (strcmp(*(atoms->atomname[i]), "NZ") == 0))
139             {
140                 atm.Cn[6] = i;
141             }
142             /* HChi flag here too */
143             else if (bHChi && (strcmp(*(atoms->atomname[i]), "NH1") == 0))
144             {
145                 atm.Cn[7] = i;
146             }
147             i++;
148         }
149
150         thisres = *(atoms->resinfo[ires].name);
151
152         /* added by grs - special case for aromatics, whose chis above 2 are
153            not real and produce rubbish output - so set back to -1 */
154         if (strcmp(thisres, "PHE") == 0 ||
155             strcmp(thisres, "TYR") == 0 ||
156             strcmp(thisres, "PTR") == 0 ||
157             strcmp(thisres, "TRP") == 0 ||
158             strcmp(thisres, "HIS") == 0 ||
159             strcmp(thisres, "HISA") == 0 ||
160             strcmp(thisres, "HISB") == 0)
161         {
162             for (ii = 5; ii <= 7; ii++)
163             {
164                 atm.Cn[ii] = -1;
165             }
166         }
167         /* end fixing aromatics */
168
169         /* Special case for Pro, has no H */
170         if (strcmp(thisres, "PRO") == 0)
171         {
172             atm.H = atm.Cn[4];
173         }
174         /* Carbon from previous residue */
175         if (prev.C != -1)
176         {
177             atm.minC = prev.C;
178         }
179         /* Alpha-carbon from previous residue */
180         if (prev.Cn[1] != -1)
181         {
182             atm.minCalpha = prev.Cn[1];
183         }
184         prev = atm;
185
186         /* Check how many dihedrals we have */
187         if ((atm.N != -1) && (atm.Cn[1] != -1) && (atm.C != -1) &&
188             (atm.O != -1) && ((atm.H != -1) || (atm.minC != -1)))
189         {
190             dl[nl].resnr     = ires+1;
191             dl[nl].atm       = atm;
192             dl[nl].atm.Cn[0] = atm.N;
193             if ((atm.Cn[3] != -1) && (atm.Cn[2] != -1) && (atm.Cn[1] != -1))
194             {
195                 nc[0]++;
196                 if (atm.Cn[4] != -1)
197                 {
198                     nc[1]++;
199                     if (atm.Cn[5] != -1)
200                     {
201                         nc[2]++;
202                         if (atm.Cn[6] != -1)
203                         {
204                             nc[3]++;
205                             if (atm.Cn[7] != -1)
206                             {
207                                 nc[4]++;
208                                 if (atm.Cn[8] != -1)
209                                 {
210                                     nc[5]++;
211                                 }
212                             }
213                         }
214                     }
215                 }
216             }
217             if ((atm.minC != -1) && (atm.minCalpha != -1))
218             {
219                 nc[6]++;
220             }
221             dl[nl].index = gmx_residuetype_get_index(rt, thisres);
222
223             sprintf(dl[nl].name, "%s%d", thisres, ires+r0);
224             nl++;
225         }
226         else if (debug)
227         {
228             fprintf(debug, "Could not find N atom but could find other atoms"
229                     " in residue %s%d\n", thisres, ires+r0);
230         }
231     }
232     fprintf(stderr, "\n");
233     fprintf(log, "\n");
234     fprintf(log, "There are %d residues with dihedrals\n", nl);
235     j = 0;
236     if (bPhi)
237     {
238         j += nl;
239     }
240     if (bPsi)
241     {
242         j += nl;
243     }
244     if (bChi)
245     {
246         for (i = 0; (i < maxchi); i++)
247         {
248             j += nc[i];
249         }
250     }
251     fprintf(log, "There are %d dihedrals\n", j);
252     fprintf(log, "Dihedral: ");
253     if (bPhi)
254     {
255         fprintf(log, " Phi  ");
256     }
257     if (bPsi)
258     {
259         fprintf(log, " Psi  ");
260     }
261     if (bChi)
262     {
263         for (i = 0; (i < maxchi); i++)
264         {
265             fprintf(log, "Chi%d  ", i+1);
266         }
267     }
268     fprintf(log, "\nNumber:   ");
269     if (bPhi)
270     {
271         fprintf(log, "%4d  ", nl);
272     }
273     if (bPsi)
274     {
275         fprintf(log, "%4d  ", nl);
276     }
277     if (bChi)
278     {
279         for (i = 0; (i < maxchi); i++)
280         {
281             fprintf(log, "%4d  ", nc[i]);
282         }
283     }
284     fprintf(log, "\n");
285
286     *nlist = nl;
287
288     return dl;
289 }
290
291 gmx_bool has_dihedral(int Dih, t_dlist *dl)
292 {
293     gmx_bool b = FALSE;
294     int      ddd;
295
296     switch (Dih)
297     {
298         case edPhi:
299             b = ((dl->atm.H != -1) && (dl->atm.N != -1) && (dl->atm.Cn[1] != -1) && (dl->atm.C != -1));
300             break;
301         case edPsi:
302             b = ((dl->atm.N != -1) && (dl->atm.Cn[1] != -1) && (dl->atm.C != -1) && (dl->atm.O != -1));
303             break;
304         case edOmega:
305             b = ((dl->atm.minCalpha != -1) && (dl->atm.minC != -1) && (dl->atm.N != -1) && (dl->atm.Cn[1] != -1));
306             break;
307         case edChi1:
308         case edChi2:
309         case edChi3:
310         case edChi4:
311         case edChi5:
312         case edChi6:
313             ddd = Dih - edChi1;
314             b   = ((dl->atm.Cn[ddd] != -1) &&  (dl->atm.Cn[ddd+1] != -1) &&
315                    (dl->atm.Cn[ddd+2] != -1) && (dl->atm.Cn[ddd+3] != -1));
316             break;
317         default:
318             pr_dlist(stdout, 1, dl, 1, 0, TRUE, TRUE, TRUE, TRUE, MAXCHI);
319             gmx_fatal(FARGS, "Non existant dihedral %d in file %s, line %d",
320                       Dih, __FILE__, __LINE__);
321     }
322     return b;
323 }
324
325 static void pr_one_ro(FILE *fp, t_dlist *dl, int nDih, real dt)
326 {
327     int k;
328     for (k = 0; k < NROT; k++)
329     {
330         fprintf(fp, "  %6.2f", dl->rot_occ[nDih][k]);
331     }
332     fprintf(fp, "\n");
333 }
334
335 static void pr_ntr_s2(FILE *fp, t_dlist *dl, int nDih, real dt)
336 {
337     fprintf(fp, "  %6.2f  %6.2f\n", (dt == 0) ? 0 : dl->ntr[nDih]/dt, dl->S2[nDih]);
338 }
339
340 void pr_dlist(FILE *fp, int nl, t_dlist dl[], real dt, int printtype,
341               gmx_bool bPhi, gmx_bool bPsi, gmx_bool bChi, gmx_bool bOmega, int maxchi)
342 {
343     int   i, Xi;
344
345     void  (*pr_props)(FILE *, t_dlist *, int, real);
346
347     /* Analysis of dihedral transitions etc */
348
349     if (printtype == edPrintST)
350     {
351         pr_props = pr_ntr_s2;
352         fprintf(stderr, "Now printing out transitions and OPs...\n");
353     }
354     else
355     {
356         pr_props = pr_one_ro;
357         fprintf(stderr, "Now printing out rotamer occupancies...\n");
358         fprintf(fp, "\nXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX\n\n");
359     }
360
361     /* change atom numbers from 0 based to 1 based */
362     for (i = 0; (i < nl); i++)
363     {
364         fprintf(fp, "Residue %s\n", dl[i].name);
365         if (printtype == edPrintST)
366         {
367             fprintf(fp, " Angle [   AI,   AJ,   AK,   AL]  #tr/ns  S^2D  \n"
368                     "--------------------------------------------\n");
369         }
370         else
371         {
372             fprintf(fp, " Angle [   AI,   AJ,   AK,   AL]  rotamers  0  g(-)  t  g(+)\n"
373                     "--------------------------------------------\n");
374         }
375         if (bPhi)
376         {
377             fprintf(fp, "   Phi [%5d,%5d,%5d,%5d]",
378                     (dl[i].atm.H == -1) ? 1+dl[i].atm.minC : 1+dl[i].atm.H,
379                     1+dl[i].atm.N, 1+dl[i].atm.Cn[1], 1+dl[i].atm.C);
380             pr_props(fp, &dl[i], edPhi, dt);
381         }
382         if (bPsi)
383         {
384             fprintf(fp, "   Psi [%5d,%5d,%5d,%5d]", 1+dl[i].atm.N, 1+dl[i].atm.Cn[1],
385                     1+dl[i].atm.C, 1+dl[i].atm.O);
386             pr_props(fp, &dl[i], edPsi, dt);
387         }
388         if (bOmega && has_dihedral(edOmega, &(dl[i])))
389         {
390             fprintf(fp, " Omega [%5d,%5d,%5d,%5d]", 1+dl[i].atm.minCalpha, 1+dl[i].atm.minC,
391                     1+dl[i].atm.N, 1+dl[i].atm.Cn[1]);
392             pr_props(fp, &dl[i], edOmega, dt);
393         }
394         for (Xi = 0; Xi < MAXCHI; Xi++)
395         {
396             if (bChi && (Xi < maxchi) && (dl[i].atm.Cn[Xi+3] != -1) )
397             {
398                 fprintf(fp, "   Chi%d[%5d,%5d,%5d,%5d]", Xi+1, 1+dl[i].atm.Cn[Xi],
399                         1+dl[i].atm.Cn[Xi+1], 1+dl[i].atm.Cn[Xi+2],
400                         1+dl[i].atm.Cn[Xi+3]);
401                 pr_props(fp, &dl[i], Xi+edChi1, dt); /* Xi+2 was wrong here */
402             }
403         }
404         fprintf(fp, "\n");
405     }
406 }
407
408
409
410 int pr_trans(FILE *fp, int nl, t_dlist dl[], real dt, int Xi)
411 {
412     /* never called at the moment */
413
414     int  i, nn, nz;
415
416     nz = 0;
417     fprintf(fp, "\\begin{table}[h]\n");
418     fprintf(fp, "\\caption{Number of dihedral transitions per nanosecond}\n");
419     fprintf(fp, "\\begin{tabular}{|l|l|}\n");
420     fprintf(fp, "\\hline\n");
421     fprintf(fp, "Residue\t&$\\chi_%d$\t\\\\\n", Xi+1);
422     for (i = 0; (i < nl); i++)
423     {
424         nn = dl[i].ntr[Xi]/dt;
425
426         if (nn == 0)
427         {
428             fprintf(fp, "%s\t&\\HL{%d}\t\\\\\n", dl[i].name, nn);
429             nz++;
430         }
431         else if (nn > 0)
432         {
433             fprintf(fp, "%s\t&\\%d\t\\\\\n", dl[i].name, nn);
434         }
435     }
436     fprintf(fp, "\\hline\n");
437     fprintf(fp, "\\end{tabular}\n");
438     fprintf(fp, "\\end{table}\n\n");
439
440     return nz;
441 }