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