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