2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
45 #include "gromacs/gmxana/gstat.h"
46 #include "gromacs/topology/residuetypes.h"
47 #include "gromacs/topology/topology.h"
48 #include "gromacs/utility/fatalerror.h"
50 std::vector<t_dlist> mk_dlist(FILE* log,
62 int nl = 0, nc[edMax];
64 // Initially, size this to all possible residues. Later it might
65 // be reduced to only handle those residues identified to contain
67 std::vector<t_dlist> dl(atoms->nres + 1);
69 prev.C = prev.Cn[1] = -1; /* Keep the compiler quiet */
70 for (i = 0; (i < edMax); i++)
77 int ires = atoms->atom[i].resind;
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++)
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))
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))
96 else if (std::strcmp(*(atoms->atomname[i]), "N") == 0)
100 else if (std::strcmp(*(atoms->atomname[i]), "C") == 0)
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))
111 else if (std::strcmp(*(atoms->atomname[i]), "CA") == 0)
115 else if (std::strcmp(*(atoms->atomname[i]), "CB") == 0)
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))
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
135 /* by grs - split the Cn[4] into 2 bits to check allowing dih to H */
137 && ((std::strcmp(*(atoms->atomname[i]), "HG") == 0)
138 || (std::strcmp(*(atoms->atomname[i]), "HG1") == 0)))
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))
149 else if ((std::strcmp(*(atoms->atomname[i]), "CZ") == 0)
150 || (std::strcmp(*(atoms->atomname[i]), "NZ") == 0))
154 /* HChi flag here too */
155 else if (bHChi && (std::strcmp(*(atoms->atomname[i]), "NH1") == 0))
162 thisres = *(atoms->resinfo[ires].name);
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)
171 for (ii = 5; ii <= 7; ii++)
176 /* end fixing aromatics */
178 /* Special case for Pro, has no H */
179 if (std::strcmp(thisres, "PRO") == 0)
183 /* Carbon from previous residue */
188 /* Alpha-carbon from previous residue */
189 if (prev.Cn[1] != -1)
191 atm.minCalpha = prev.Cn[1];
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)))
199 dl[nl].resnr = ires + 1;
201 dl[nl].atm.Cn[0] = atm.N;
202 if ((atm.Cn[3] != -1) && (atm.Cn[2] != -1) && (atm.Cn[1] != -1))
226 if ((atm.minC != -1) && (atm.minCalpha != -1))
230 dl[nl].index = rt->indexFromResidueName(thisres);
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)
238 "Unknown residue %s when searching for residue type.\n"
239 "Maybe you need to add a custom residue in residuetypes.dat.",
243 sprintf(dl[nl].name, "%s%d", thisres, ires + r0);
249 "Could not find N atom but could find other atoms"
250 " in residue %s%d\n",
255 // Leave only the residues that were recognized to contain dihedrals
258 fprintf(stderr, "\n");
260 fprintf(log, "There are %d residues with dihedrals\n", nl);
272 for (i = 0; (i < maxchi); i++)
277 fprintf(log, "There are %d dihedrals\n", j);
278 fprintf(log, "Dihedral: ");
281 fprintf(log, " Phi ");
285 fprintf(log, " Psi ");
289 for (i = 0; (i < maxchi); i++)
291 fprintf(log, "Chi%d ", i + 1);
294 fprintf(log, "\nNumber: ");
297 fprintf(log, "%4d ", nl);
301 fprintf(log, "%4d ", nl);
305 for (i = 0; (i < maxchi); i++)
307 fprintf(log, "%4d ", nc[i]);
315 gmx_bool has_dihedral(int Dih, const t_dlist& dl)
323 b = ((dl.atm.H != -1) && (dl.atm.N != -1) && (dl.atm.Cn[1] != -1) && (dl.atm.C != -1));
326 b = ((dl.atm.N != -1) && (dl.atm.Cn[1] != -1) && (dl.atm.C != -1) && (dl.atm.O != -1));
329 b = ((dl.atm.minCalpha != -1) && (dl.atm.minC != -1) && (dl.atm.N != -1)
330 && (dl.atm.Cn[1] != -1));
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));
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__);
349 static void pr_one_ro(FILE* fp, const t_dlist& dl, int nDih, real gmx_unused dt)
352 for (k = 0; k < NROT; k++)
354 fprintf(fp, " %6.2f", dl.rot_occ[nDih][k]);
359 static void pr_ntr_s2(FILE* fp, const t_dlist& dl, int nDih, real dt)
361 fprintf(fp, " %6.2f %6.2f\n", (dt == 0) ? 0 : dl.ntr[nDih] / dt, dl.S2[nDih]);
364 void pr_dlist(FILE* fp,
365 gmx::ArrayRef<const t_dlist> dlist,
374 void (*pr_props)(FILE*, const t_dlist&, int, real);
376 /* Analysis of dihedral transitions etc */
378 if (printtype == edPrintST)
380 pr_props = pr_ntr_s2;
381 fprintf(stderr, "Now printing out transitions and OPs...\n");
385 pr_props = pr_one_ro;
386 fprintf(stderr, "Now printing out rotamer occupancies...\n");
387 fprintf(fp, "\nXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX\n\n");
390 /* change atom numbers from 0 based to 1 based */
391 for (const auto& dihedral : dlist)
393 fprintf(fp, "Residue %s\n", dihedral.name);
394 if (printtype == edPrintST)
397 " Angle [ AI, AJ, AK, AL] #tr/ns S^2D \n"
398 "--------------------------------------------\n");
403 " Angle [ AI, AJ, AK, AL] rotamers 0 g(-) t g(+)\n"
404 "--------------------------------------------\n");
409 " Phi [%5d,%5d,%5d,%5d]",
410 (dihedral.atm.H == -1) ? 1 + dihedral.atm.minC : 1 + dihedral.atm.H,
412 1 + dihedral.atm.Cn[1],
414 pr_props(fp, dihedral, edPhi, dt);
419 " Psi [%5d,%5d,%5d,%5d]",
421 1 + dihedral.atm.Cn[1],
424 pr_props(fp, dihedral, edPsi, dt);
426 if (bOmega && has_dihedral(edOmega, dihedral))
429 " Omega [%5d,%5d,%5d,%5d]",
430 1 + dihedral.atm.minCalpha,
431 1 + dihedral.atm.minC,
433 1 + dihedral.atm.Cn[1]);
434 pr_props(fp, dihedral, edOmega, dt);
436 for (int Xi = 0; Xi < MAXCHI; Xi++)
438 if (bChi && (Xi < maxchi) && (dihedral.atm.Cn[Xi + 3] != -1))
441 " Chi%d[%5d,%5d,%5d,%5d]",
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 */