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