6d9fbacc8e3b85e4e45aedd0f684fa7426783a46
[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,2021, 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 <vector>
44
45 #include "gromacs/gmxana/gstat.h"
46 #include "gromacs/topology/residuetypes.h"
47 #include "gromacs/topology/topology.h"
48 #include "gromacs/utility/fatalerror.h"
49
50 std::vector<t_dlist> mk_dlist(FILE*          log,
51                               const t_atoms* atoms,
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     // Initially, size this to all possible residues. Later it might
65     // be reduced to only handle those residues identified to contain
66     // dihedrals.
67     std::vector<t_dlist> dl(atoms->nres + 1);
68
69     prev.C = prev.Cn[1] = -1; /* Keep the compiler quiet */
70     for (i = 0; (i < edMax); i++)
71     {
72         nc[i] = 0;
73     }
74     i = 0;
75     while (i < atoms->nr)
76     {
77         int ires = atoms->atom[i].resind;
78
79         /* Initiate all atom numbers to -1 */
80         atm.minC = atm.H = atm.N = atm.C = atm.O = atm.minCalpha = -1;
81         for (j = 0; (j < MAXCHI + 3); j++)
82         {
83             atm.Cn[j] = -1;
84         }
85
86         /* Look for atoms in this residue */
87         /* maybe should allow for chis to hydrogens? */
88         while ((i < atoms->nr) && (atoms->atom[i].resind == ires))
89         {
90             if ((std::strcmp(*(atoms->atomname[i]), "H") == 0)
91                 || (std::strcmp(*(atoms->atomname[i]), "H1") == 0)
92                 || (std::strcmp(*(atoms->atomname[i]), "HN") == 0))
93             {
94                 atm.H = i;
95             }
96             else if (std::strcmp(*(atoms->atomname[i]), "N") == 0)
97             {
98                 atm.N = i;
99             }
100             else if (std::strcmp(*(atoms->atomname[i]), "C") == 0)
101             {
102                 atm.C = i;
103             }
104             else if ((std::strcmp(*(atoms->atomname[i]), "O") == 0)
105                      || (std::strcmp(*(atoms->atomname[i]), "O1") == 0)
106                      || (std::strcmp(*(atoms->atomname[i]), "OC1") == 0)
107                      || (std::strcmp(*(atoms->atomname[i]), "OT1") == 0))
108             {
109                 atm.O = i;
110             }
111             else if (std::strcmp(*(atoms->atomname[i]), "CA") == 0)
112             {
113                 atm.Cn[1] = i;
114             }
115             else if (std::strcmp(*(atoms->atomname[i]), "CB") == 0)
116             {
117                 atm.Cn[2] = i;
118             }
119             else if ((std::strcmp(*(atoms->atomname[i]), "CG") == 0)
120                      || (std::strcmp(*(atoms->atomname[i]), "CG1") == 0)
121                      || (std::strcmp(*(atoms->atomname[i]), "OG") == 0)
122                      || (std::strcmp(*(atoms->atomname[i]), "OG1") == 0)
123                      || (std::strcmp(*(atoms->atomname[i]), "SG") == 0))
124             {
125                 atm.Cn[3] = i;
126             }
127             else if ((std::strcmp(*(atoms->atomname[i]), "CD") == 0)
128                      || (std::strcmp(*(atoms->atomname[i]), "CD1") == 0)
129                      || (std::strcmp(*(atoms->atomname[i]), "SD") == 0)
130                      || (std::strcmp(*(atoms->atomname[i]), "OD1") == 0)
131                      || (std::strcmp(*(atoms->atomname[i]), "ND1") == 0))
132             { // NOLINT bugprone-branch-clone
133                 atm.Cn[4] = i;
134             }
135             /* by grs - split the Cn[4] into 2 bits to check allowing dih to H */
136             else if (bHChi
137                      && ((std::strcmp(*(atoms->atomname[i]), "HG") == 0)
138                          || (std::strcmp(*(atoms->atomname[i]), "HG1") == 0)))
139             {
140                 atm.Cn[4] = i;
141             }
142             else if ((std::strcmp(*(atoms->atomname[i]), "CE") == 0)
143                      || (std::strcmp(*(atoms->atomname[i]), "CE1") == 0)
144                      || (std::strcmp(*(atoms->atomname[i]), "OE1") == 0)
145                      || (std::strcmp(*(atoms->atomname[i]), "NE") == 0))
146             {
147                 atm.Cn[5] = i;
148             }
149             else if ((std::strcmp(*(atoms->atomname[i]), "CZ") == 0)
150                      || (std::strcmp(*(atoms->atomname[i]), "NZ") == 0))
151             {
152                 atm.Cn[6] = i;
153             }
154             /* HChi flag here too */
155             else if (bHChi && (std::strcmp(*(atoms->atomname[i]), "NH1") == 0))
156             {
157                 atm.Cn[7] = i;
158             }
159             i++;
160         }
161
162         thisres = *(atoms->resinfo[ires].name);
163
164         /* added by grs - special case for aromatics, whose chis above 2 are
165            not real and produce rubbish output - so set back to -1 */
166         if (std::strcmp(thisres, "PHE") == 0 || std::strcmp(thisres, "TYR") == 0
167             || std::strcmp(thisres, "PTR") == 0 || std::strcmp(thisres, "TRP") == 0
168             || std::strcmp(thisres, "HIS") == 0 || std::strcmp(thisres, "HISA") == 0
169             || std::strcmp(thisres, "HISB") == 0)
170         {
171             for (ii = 5; ii <= 7; ii++)
172             {
173                 atm.Cn[ii] = -1;
174             }
175         }
176         /* end fixing aromatics */
177
178         /* Special case for Pro, has no H */
179         if (std::strcmp(thisres, "PRO") == 0)
180         {
181             atm.H = atm.Cn[4];
182         }
183         /* Carbon from previous residue */
184         if (prev.C != -1)
185         {
186             atm.minC = prev.C;
187         }
188         /* Alpha-carbon from previous residue */
189         if (prev.Cn[1] != -1)
190         {
191             atm.minCalpha = prev.Cn[1];
192         }
193         prev = atm;
194
195         /* Check how many dihedrals we have */
196         if ((atm.N != -1) && (atm.Cn[1] != -1) && (atm.C != -1) && (atm.O != -1)
197             && ((atm.H != -1) || (atm.minC != -1)))
198         {
199             dl[nl].resnr     = ires + 1;
200             dl[nl].atm       = atm;
201             dl[nl].atm.Cn[0] = atm.N;
202             if ((atm.Cn[3] != -1) && (atm.Cn[2] != -1) && (atm.Cn[1] != -1))
203             {
204                 nc[0]++;
205                 if (atm.Cn[4] != -1)
206                 {
207                     nc[1]++;
208                     if (atm.Cn[5] != -1)
209                     {
210                         nc[2]++;
211                         if (atm.Cn[6] != -1)
212                         {
213                             nc[3]++;
214                             if (atm.Cn[7] != -1)
215                             {
216                                 nc[4]++;
217                                 if (atm.Cn[8] != -1)
218                                 {
219                                     nc[5]++;
220                                 }
221                             }
222                         }
223                     }
224                 }
225             }
226             if ((atm.minC != -1) && (atm.minCalpha != -1))
227             {
228                 nc[6]++;
229             }
230             dl[nl].index = rt->indexFromResidueName(thisres);
231
232             /* Prevent use of unknown residues. If one adds a custom residue to
233              * residuetypes.dat but somehow loses it, changes it, or does analysis on
234              * another machine, the residue type will be unknown. */
235             if (dl[nl].index == -1)
236             {
237                 gmx_fatal(FARGS,
238                           "Unknown residue %s when searching for residue type.\n"
239                           "Maybe you need to add a custom residue in residuetypes.dat.",
240                           thisres);
241             }
242
243             sprintf(dl[nl].name, "%s%d", thisres, ires + r0);
244             nl++;
245         }
246         else if (debug)
247         {
248             fprintf(debug,
249                     "Could not find N atom but could find other atoms"
250                     " in residue %s%d\n",
251                     thisres,
252                     ires + r0);
253         }
254     }
255     // Leave only the residues that were recognized to contain dihedrals
256     dl.resize(nl);
257
258     fprintf(stderr, "\n");
259     fprintf(log, "\n");
260     fprintf(log, "There are %d residues with dihedrals\n", nl);
261     j = 0;
262     if (bPhi)
263     {
264         j += nl;
265     }
266     if (bPsi)
267     {
268         j += nl;
269     }
270     if (bChi)
271     {
272         for (i = 0; (i < maxchi); i++)
273         {
274             j += nc[i];
275         }
276     }
277     fprintf(log, "There are %d dihedrals\n", j);
278     fprintf(log, "Dihedral: ");
279     if (bPhi)
280     {
281         fprintf(log, " Phi  ");
282     }
283     if (bPsi)
284     {
285         fprintf(log, " Psi  ");
286     }
287     if (bChi)
288     {
289         for (i = 0; (i < maxchi); i++)
290         {
291             fprintf(log, "Chi%d  ", i + 1);
292         }
293     }
294     fprintf(log, "\nNumber:   ");
295     if (bPhi)
296     {
297         fprintf(log, "%4d  ", nl);
298     }
299     if (bPsi)
300     {
301         fprintf(log, "%4d  ", nl);
302     }
303     if (bChi)
304     {
305         for (i = 0; (i < maxchi); i++)
306         {
307             fprintf(log, "%4d  ", nc[i]);
308         }
309     }
310     fprintf(log, "\n");
311
312     return dl;
313 }
314
315 gmx_bool has_dihedral(int Dih, const t_dlist& dl)
316 {
317     gmx_bool b = FALSE;
318     int      ddd;
319
320     switch (Dih)
321     {
322         case edPhi:
323             b = ((dl.atm.H != -1) && (dl.atm.N != -1) && (dl.atm.Cn[1] != -1) && (dl.atm.C != -1));
324             break;
325         case edPsi:
326             b = ((dl.atm.N != -1) && (dl.atm.Cn[1] != -1) && (dl.atm.C != -1) && (dl.atm.O != -1));
327             break;
328         case edOmega:
329             b = ((dl.atm.minCalpha != -1) && (dl.atm.minC != -1) && (dl.atm.N != -1)
330                  && (dl.atm.Cn[1] != -1));
331             break;
332         case edChi1:
333         case edChi2:
334         case edChi3:
335         case edChi4:
336         case edChi5:
337         case edChi6:
338             ddd = Dih - edChi1;
339             b = ((dl.atm.Cn[ddd] != -1) && (dl.atm.Cn[ddd + 1] != -1) && (dl.atm.Cn[ddd + 2] != -1)
340                  && (dl.atm.Cn[ddd + 3] != -1));
341             break;
342         default:
343             pr_dlist(stdout, gmx::constArrayRefFromArray(&dl, 1), 1, 0, TRUE, TRUE, TRUE, TRUE, MAXCHI);
344             gmx_fatal(FARGS, "Non existent dihedral %d in file %s, line %d", Dih, __FILE__, __LINE__);
345     }
346     return b;
347 }
348
349 static void pr_one_ro(FILE* fp, const t_dlist& dl, int nDih, real gmx_unused dt)
350 {
351     int k;
352     for (k = 0; k < NROT; k++)
353     {
354         fprintf(fp, "  %6.2f", dl.rot_occ[nDih][k]);
355     }
356     fprintf(fp, "\n");
357 }
358
359 static void pr_ntr_s2(FILE* fp, const t_dlist& dl, int nDih, real dt)
360 {
361     fprintf(fp, "  %6.2f  %6.2f\n", (dt == 0) ? 0 : dl.ntr[nDih] / dt, dl.S2[nDih]);
362 }
363
364 void pr_dlist(FILE*                        fp,
365               gmx::ArrayRef<const t_dlist> dlist,
366               real                         dt,
367               int                          printtype,
368               gmx_bool                     bPhi,
369               gmx_bool                     bPsi,
370               gmx_bool                     bChi,
371               gmx_bool                     bOmega,
372               int                          maxchi)
373 {
374     void (*pr_props)(FILE*, const 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 (const auto& dihedral : dlist)
392     {
393         fprintf(fp, "Residue %s\n", dihedral.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                     (dihedral.atm.H == -1) ? 1 + dihedral.atm.minC : 1 + dihedral.atm.H,
411                     1 + dihedral.atm.N,
412                     1 + dihedral.atm.Cn[1],
413                     1 + dihedral.atm.C);
414             pr_props(fp, dihedral, edPhi, dt);
415         }
416         if (bPsi)
417         {
418             fprintf(fp,
419                     "   Psi [%5d,%5d,%5d,%5d]",
420                     1 + dihedral.atm.N,
421                     1 + dihedral.atm.Cn[1],
422                     1 + dihedral.atm.C,
423                     1 + dihedral.atm.O);
424             pr_props(fp, dihedral, edPsi, dt);
425         }
426         if (bOmega && has_dihedral(edOmega, dihedral))
427         {
428             fprintf(fp,
429                     " Omega [%5d,%5d,%5d,%5d]",
430                     1 + dihedral.atm.minCalpha,
431                     1 + dihedral.atm.minC,
432                     1 + dihedral.atm.N,
433                     1 + dihedral.atm.Cn[1]);
434             pr_props(fp, dihedral, edOmega, dt);
435         }
436         for (int Xi = 0; Xi < MAXCHI; Xi++)
437         {
438             if (bChi && (Xi < maxchi) && (dihedral.atm.Cn[Xi + 3] != -1))
439             {
440                 fprintf(fp,
441                         "   Chi%d[%5d,%5d,%5d,%5d]",
442                         Xi + 1,
443                         1 + dihedral.atm.Cn[Xi],
444                         1 + dihedral.atm.Cn[Xi + 1],
445                         1 + dihedral.atm.Cn[Xi + 2],
446                         1 + dihedral.atm.Cn[Xi + 3]);
447                 pr_props(fp, dihedral, Xi + edChi1, dt); /* Xi+2 was wrong here */
448             }
449         }
450         fprintf(fp, "\n");
451     }
452 }