Apply re-formatting to C++ in src/ tree.
[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,2016,2018 by the GROMACS development team.
7  * Copyright (c) 2019,2020, by the GROMACS development team, led by
8  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9  * and including many others, as listed in the AUTHORS file in the
10  * top-level source directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #include "gmxpre.h"
39
40 #include <cstdlib>
41 #include <cstring>
42
43 #include "gromacs/gmxana/gstat.h"
44 #include "gromacs/topology/residuetypes.h"
45 #include "gromacs/topology/topology.h"
46 #include "gromacs/utility/fatalerror.h"
47 #include "gromacs/utility/smalloc.h"
48
49 t_dlist* mk_dlist(FILE*          log,
50                   const t_atoms* atoms,
51                   int*           nlist,
52                   gmx_bool       bPhi,
53                   gmx_bool       bPsi,
54                   gmx_bool       bChi,
55                   gmx_bool       bHChi,
56                   int            maxchi,
57                   int            r0,
58                   ResidueType*   rt)
59 {
60     int       i, j, ii;
61     t_dihatms atm, prev;
62     int       nl = 0, nc[edMax];
63     char*     thisres;
64     t_dlist*  dl;
65
66     snew(dl, atoms->nres + 1);
67     prev.C = prev.Cn[1] = -1; /* Keep the compiler quiet */
68     for (i = 0; (i < edMax); i++)
69     {
70         nc[i] = 0;
71     }
72     i = 0;
73     while (i < atoms->nr)
74     {
75         int ires = atoms->atom[i].resind;
76
77         /* Initiate all atom numbers to -1 */
78         atm.minC = atm.H = atm.N = atm.C = atm.O = atm.minCalpha = -1;
79         for (j = 0; (j < MAXCHI + 3); j++)
80         {
81             atm.Cn[j] = -1;
82         }
83
84         /* Look for atoms in this residue */
85         /* maybe should allow for chis to hydrogens? */
86         while ((i < atoms->nr) && (atoms->atom[i].resind == ires))
87         {
88             if ((std::strcmp(*(atoms->atomname[i]), "H") == 0)
89                 || (std::strcmp(*(atoms->atomname[i]), "H1") == 0)
90                 || (std::strcmp(*(atoms->atomname[i]), "HN") == 0))
91             {
92                 atm.H = i;
93             }
94             else if (std::strcmp(*(atoms->atomname[i]), "N") == 0)
95             {
96                 atm.N = i;
97             }
98             else if (std::strcmp(*(atoms->atomname[i]), "C") == 0)
99             {
100                 atm.C = i;
101             }
102             else if ((std::strcmp(*(atoms->atomname[i]), "O") == 0)
103                      || (std::strcmp(*(atoms->atomname[i]), "O1") == 0)
104                      || (std::strcmp(*(atoms->atomname[i]), "OC1") == 0)
105                      || (std::strcmp(*(atoms->atomname[i]), "OT1") == 0))
106             {
107                 atm.O = i;
108             }
109             else if (std::strcmp(*(atoms->atomname[i]), "CA") == 0)
110             {
111                 atm.Cn[1] = i;
112             }
113             else if (std::strcmp(*(atoms->atomname[i]), "CB") == 0)
114             {
115                 atm.Cn[2] = i;
116             }
117             else if ((std::strcmp(*(atoms->atomname[i]), "CG") == 0)
118                      || (std::strcmp(*(atoms->atomname[i]), "CG1") == 0)
119                      || (std::strcmp(*(atoms->atomname[i]), "OG") == 0)
120                      || (std::strcmp(*(atoms->atomname[i]), "OG1") == 0)
121                      || (std::strcmp(*(atoms->atomname[i]), "SG") == 0))
122             {
123                 atm.Cn[3] = i;
124             }
125             else if ((std::strcmp(*(atoms->atomname[i]), "CD") == 0)
126                      || (std::strcmp(*(atoms->atomname[i]), "CD1") == 0)
127                      || (std::strcmp(*(atoms->atomname[i]), "SD") == 0)
128                      || (std::strcmp(*(atoms->atomname[i]), "OD1") == 0)
129                      || (std::strcmp(*(atoms->atomname[i]), "ND1") == 0))
130             { // NOLINT bugprone-branch-clone
131                 atm.Cn[4] = i;
132             }
133             /* by grs - split the Cn[4] into 2 bits to check allowing dih to H */
134             else if (bHChi
135                      && ((std::strcmp(*(atoms->atomname[i]), "HG") == 0)
136                          || (std::strcmp(*(atoms->atomname[i]), "HG1") == 0)))
137             {
138                 atm.Cn[4] = i;
139             }
140             else if ((std::strcmp(*(atoms->atomname[i]), "CE") == 0)
141                      || (std::strcmp(*(atoms->atomname[i]), "CE1") == 0)
142                      || (std::strcmp(*(atoms->atomname[i]), "OE1") == 0)
143                      || (std::strcmp(*(atoms->atomname[i]), "NE") == 0))
144             {
145                 atm.Cn[5] = i;
146             }
147             else if ((std::strcmp(*(atoms->atomname[i]), "CZ") == 0)
148                      || (std::strcmp(*(atoms->atomname[i]), "NZ") == 0))
149             {
150                 atm.Cn[6] = i;
151             }
152             /* HChi flag here too */
153             else if (bHChi && (std::strcmp(*(atoms->atomname[i]), "NH1") == 0))
154             {
155                 atm.Cn[7] = i;
156             }
157             i++;
158         }
159
160         thisres = *(atoms->resinfo[ires].name);
161
162         /* added by grs - special case for aromatics, whose chis above 2 are
163            not real and produce rubbish output - so set back to -1 */
164         if (std::strcmp(thisres, "PHE") == 0 || std::strcmp(thisres, "TYR") == 0
165             || std::strcmp(thisres, "PTR") == 0 || std::strcmp(thisres, "TRP") == 0
166             || std::strcmp(thisres, "HIS") == 0 || std::strcmp(thisres, "HISA") == 0
167             || std::strcmp(thisres, "HISB") == 0)
168         {
169             for (ii = 5; ii <= 7; ii++)
170             {
171                 atm.Cn[ii] = -1;
172             }
173         }
174         /* end fixing aromatics */
175
176         /* Special case for Pro, has no H */
177         if (std::strcmp(thisres, "PRO") == 0)
178         {
179             atm.H = atm.Cn[4];
180         }
181         /* Carbon from previous residue */
182         if (prev.C != -1)
183         {
184             atm.minC = prev.C;
185         }
186         /* Alpha-carbon from previous residue */
187         if (prev.Cn[1] != -1)
188         {
189             atm.minCalpha = prev.Cn[1];
190         }
191         prev = atm;
192
193         /* Check how many dihedrals we have */
194         if ((atm.N != -1) && (atm.Cn[1] != -1) && (atm.C != -1) && (atm.O != -1)
195             && ((atm.H != -1) || (atm.minC != -1)))
196         {
197             dl[nl].resnr     = ires + 1;
198             dl[nl].atm       = atm;
199             dl[nl].atm.Cn[0] = atm.N;
200             if ((atm.Cn[3] != -1) && (atm.Cn[2] != -1) && (atm.Cn[1] != -1))
201             {
202                 nc[0]++;
203                 if (atm.Cn[4] != -1)
204                 {
205                     nc[1]++;
206                     if (atm.Cn[5] != -1)
207                     {
208                         nc[2]++;
209                         if (atm.Cn[6] != -1)
210                         {
211                             nc[3]++;
212                             if (atm.Cn[7] != -1)
213                             {
214                                 nc[4]++;
215                                 if (atm.Cn[8] != -1)
216                                 {
217                                     nc[5]++;
218                                 }
219                             }
220                         }
221                     }
222                 }
223             }
224             if ((atm.minC != -1) && (atm.minCalpha != -1))
225             {
226                 nc[6]++;
227             }
228             dl[nl].index = rt->indexFromResidueName(thisres);
229
230             /* Prevent seg fault from unknown residues. If one adds a custom residue to
231              * residuetypes.dat but somehow loses it, changes it, or does analysis on
232              * another machine, the residue type will be unknown. */
233             if (dl[nl].index == -1)
234             {
235                 gmx_fatal(FARGS,
236                           "Unknown residue %s when searching for residue type.\n"
237                           "Maybe you need to add a custom residue in residuetypes.dat.",
238                           thisres);
239             }
240
241             sprintf(dl[nl].name, "%s%d", thisres, ires + r0);
242             nl++;
243         }
244         else if (debug)
245         {
246             fprintf(debug,
247                     "Could not find N atom but could find other atoms"
248                     " in residue %s%d\n",
249                     thisres,
250                     ires + r0);
251         }
252     }
253     fprintf(stderr, "\n");
254     fprintf(log, "\n");
255     fprintf(log, "There are %d residues with dihedrals\n", nl);
256     j = 0;
257     if (bPhi)
258     {
259         j += nl;
260     }
261     if (bPsi)
262     {
263         j += nl;
264     }
265     if (bChi)
266     {
267         for (i = 0; (i < maxchi); i++)
268         {
269             j += nc[i];
270         }
271     }
272     fprintf(log, "There are %d dihedrals\n", j);
273     fprintf(log, "Dihedral: ");
274     if (bPhi)
275     {
276         fprintf(log, " Phi  ");
277     }
278     if (bPsi)
279     {
280         fprintf(log, " Psi  ");
281     }
282     if (bChi)
283     {
284         for (i = 0; (i < maxchi); i++)
285         {
286             fprintf(log, "Chi%d  ", i + 1);
287         }
288     }
289     fprintf(log, "\nNumber:   ");
290     if (bPhi)
291     {
292         fprintf(log, "%4d  ", nl);
293     }
294     if (bPsi)
295     {
296         fprintf(log, "%4d  ", nl);
297     }
298     if (bChi)
299     {
300         for (i = 0; (i < maxchi); i++)
301         {
302             fprintf(log, "%4d  ", nc[i]);
303         }
304     }
305     fprintf(log, "\n");
306
307     *nlist = nl;
308
309     return dl;
310 }
311
312 gmx_bool has_dihedral(int Dih, t_dlist* dl)
313 {
314     gmx_bool b = FALSE;
315     int      ddd;
316
317     switch (Dih)
318     {
319         case edPhi:
320             b = ((dl->atm.H != -1) && (dl->atm.N != -1) && (dl->atm.Cn[1] != -1) && (dl->atm.C != -1));
321             break;
322         case edPsi:
323             b = ((dl->atm.N != -1) && (dl->atm.Cn[1] != -1) && (dl->atm.C != -1) && (dl->atm.O != -1));
324             break;
325         case edOmega:
326             b = ((dl->atm.minCalpha != -1) && (dl->atm.minC != -1) && (dl->atm.N != -1)
327                  && (dl->atm.Cn[1] != -1));
328             break;
329         case edChi1:
330         case edChi2:
331         case edChi3:
332         case edChi4:
333         case edChi5:
334         case edChi6:
335             ddd = Dih - edChi1;
336             b   = ((dl->atm.Cn[ddd] != -1) && (dl->atm.Cn[ddd + 1] != -1)
337                  && (dl->atm.Cn[ddd + 2] != -1) && (dl->atm.Cn[ddd + 3] != -1));
338             break;
339         default:
340             pr_dlist(stdout, 1, dl, 1, 0, TRUE, TRUE, TRUE, TRUE, MAXCHI);
341             gmx_fatal(FARGS, "Non existent dihedral %d in file %s, line %d", Dih, __FILE__, __LINE__);
342     }
343     return b;
344 }
345
346 static void pr_one_ro(FILE* fp, t_dlist* dl, int nDih, real gmx_unused dt)
347 {
348     int k;
349     for (k = 0; k < NROT; k++)
350     {
351         fprintf(fp, "  %6.2f", dl->rot_occ[nDih][k]);
352     }
353     fprintf(fp, "\n");
354 }
355
356 static void pr_ntr_s2(FILE* fp, t_dlist* dl, int nDih, real dt)
357 {
358     fprintf(fp, "  %6.2f  %6.2f\n", (dt == 0) ? 0 : dl->ntr[nDih] / dt, dl->S2[nDih]);
359 }
360
361 void pr_dlist(FILE*    fp,
362               int      nl,
363               t_dlist  dl[],
364               real     dt,
365               int      printtype,
366               gmx_bool bPhi,
367               gmx_bool bPsi,
368               gmx_bool bChi,
369               gmx_bool bOmega,
370               int      maxchi)
371 {
372     int i, Xi;
373
374     void (*pr_props)(FILE*, t_dlist*, int, real);
375
376     /* Analysis of dihedral transitions etc */
377
378     if (printtype == edPrintST)
379     {
380         pr_props = pr_ntr_s2;
381         fprintf(stderr, "Now printing out transitions and OPs...\n");
382     }
383     else
384     {
385         pr_props = pr_one_ro;
386         fprintf(stderr, "Now printing out rotamer occupancies...\n");
387         fprintf(fp, "\nXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX\n\n");
388     }
389
390     /* change atom numbers from 0 based to 1 based */
391     for (i = 0; (i < nl); i++)
392     {
393         fprintf(fp, "Residue %s\n", dl[i].name);
394         if (printtype == edPrintST)
395         {
396             fprintf(fp,
397                     " Angle [   AI,   AJ,   AK,   AL]  #tr/ns  S^2D  \n"
398                     "--------------------------------------------\n");
399         }
400         else
401         {
402             fprintf(fp,
403                     " Angle [   AI,   AJ,   AK,   AL]  rotamers  0  g(-)  t  g(+)\n"
404                     "--------------------------------------------\n");
405         }
406         if (bPhi)
407         {
408             fprintf(fp,
409                     "   Phi [%5d,%5d,%5d,%5d]",
410                     (dl[i].atm.H == -1) ? 1 + dl[i].atm.minC : 1 + dl[i].atm.H,
411                     1 + dl[i].atm.N,
412                     1 + dl[i].atm.Cn[1],
413                     1 + dl[i].atm.C);
414             pr_props(fp, &dl[i], edPhi, dt);
415         }
416         if (bPsi)
417         {
418             fprintf(fp,
419                     "   Psi [%5d,%5d,%5d,%5d]",
420                     1 + dl[i].atm.N,
421                     1 + dl[i].atm.Cn[1],
422                     1 + dl[i].atm.C,
423                     1 + dl[i].atm.O);
424             pr_props(fp, &dl[i], edPsi, dt);
425         }
426         if (bOmega && has_dihedral(edOmega, &(dl[i])))
427         {
428             fprintf(fp,
429                     " Omega [%5d,%5d,%5d,%5d]",
430                     1 + dl[i].atm.minCalpha,
431                     1 + dl[i].atm.minC,
432                     1 + dl[i].atm.N,
433                     1 + dl[i].atm.Cn[1]);
434             pr_props(fp, &dl[i], edOmega, dt);
435         }
436         for (Xi = 0; Xi < MAXCHI; Xi++)
437         {
438             if (bChi && (Xi < maxchi) && (dl[i].atm.Cn[Xi + 3] != -1))
439             {
440                 fprintf(fp,
441                         "   Chi%d[%5d,%5d,%5d,%5d]",
442                         Xi + 1,
443                         1 + dl[i].atm.Cn[Xi],
444                         1 + dl[i].atm.Cn[Xi + 1],
445                         1 + dl[i].atm.Cn[Xi + 2],
446                         1 + dl[i].atm.Cn[Xi + 3]);
447                 pr_props(fp, &dl[i], Xi + edChi1, dt); /* Xi+2 was wrong here */
448             }
449         }
450         fprintf(fp, "\n");
451     }
452 }
453
454
455 int pr_trans(FILE* fp, int nl, t_dlist dl[], real dt, int Xi)
456 {
457     /* never called at the moment */
458
459     int i, nn, nz;
460
461     nz = 0;
462     fprintf(fp, "\\begin{table}[h]\n");
463     fprintf(fp, "\\caption{Number of dihedral transitions per nanosecond}\n");
464     fprintf(fp, "\\begin{tabular}{|l|l|}\n");
465     fprintf(fp, "\\hline\n");
466     fprintf(fp, "Residue\t&$\\chi_%d$\t\\\\\n", Xi + 1);
467     for (i = 0; (i < nl); i++)
468     {
469         nn = static_cast<int>(dl[i].ntr[Xi] / dt);
470
471         if (nn == 0)
472         {
473             fprintf(fp, "%s\t&\\HL{%d}\t\\\\\n", dl[i].name, nn);
474             nz++;
475         }
476         else if (nn > 0)
477         {
478             fprintf(fp, "%s\t&\\%d\t\\\\\n", dl[i].name, nn);
479         }
480     }
481     fprintf(fp, "\\hline\n");
482     fprintf(fp, "\\end{tabular}\n");
483     fprintf(fp, "\\end{table}\n\n");
484
485     return nz;
486 }