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