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