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/topology.h"
47 #include "gromacs/utility/fatalerror.h"
49 std::vector<t_dlist> mk_dlist(FILE* log,
57 const ResidueTypeMap& residueTypeMap)
61 int nl = 0, nc[edMax];
63 // Initially, size this to all possible residues. Later it might
64 // be reduced to only handle those residues identified to contain
66 std::vector<t_dlist> dl(atoms->nres + 1);
68 prev.C = prev.Cn[1] = -1; /* Keep the compiler quiet */
69 for (i = 0; (i < edMax); i++)
76 int ires = atoms->atom[i].resind;
78 /* Initiate all atom numbers to -1 */
79 atm.minC = atm.H = atm.N = atm.C = atm.O = atm.minCalpha = -1;
80 for (j = 0; (j < MAXCHI + 3); j++)
85 /* Look for atoms in this residue */
86 /* maybe should allow for chis to hydrogens? */
87 while ((i < atoms->nr) && (atoms->atom[i].resind == ires))
89 if ((std::strcmp(*(atoms->atomname[i]), "H") == 0)
90 || (std::strcmp(*(atoms->atomname[i]), "H1") == 0)
91 || (std::strcmp(*(atoms->atomname[i]), "HN") == 0))
95 else if (std::strcmp(*(atoms->atomname[i]), "N") == 0)
99 else if (std::strcmp(*(atoms->atomname[i]), "C") == 0)
103 else if ((std::strcmp(*(atoms->atomname[i]), "O") == 0)
104 || (std::strcmp(*(atoms->atomname[i]), "O1") == 0)
105 || (std::strcmp(*(atoms->atomname[i]), "OC1") == 0)
106 || (std::strcmp(*(atoms->atomname[i]), "OT1") == 0))
110 else if (std::strcmp(*(atoms->atomname[i]), "CA") == 0)
114 else if (std::strcmp(*(atoms->atomname[i]), "CB") == 0)
118 else if ((std::strcmp(*(atoms->atomname[i]), "CG") == 0)
119 || (std::strcmp(*(atoms->atomname[i]), "CG1") == 0)
120 || (std::strcmp(*(atoms->atomname[i]), "OG") == 0)
121 || (std::strcmp(*(atoms->atomname[i]), "OG1") == 0)
122 || (std::strcmp(*(atoms->atomname[i]), "SG") == 0))
126 else if ((std::strcmp(*(atoms->atomname[i]), "CD") == 0)
127 || (std::strcmp(*(atoms->atomname[i]), "CD1") == 0)
128 || (std::strcmp(*(atoms->atomname[i]), "SD") == 0)
129 || (std::strcmp(*(atoms->atomname[i]), "OD1") == 0)
130 || (std::strcmp(*(atoms->atomname[i]), "ND1") == 0))
131 { // NOLINT bugprone-branch-clone
134 /* by grs - split the Cn[4] into 2 bits to check allowing dih to H */
136 && ((std::strcmp(*(atoms->atomname[i]), "HG") == 0)
137 || (std::strcmp(*(atoms->atomname[i]), "HG1") == 0)))
141 else if ((std::strcmp(*(atoms->atomname[i]), "CE") == 0)
142 || (std::strcmp(*(atoms->atomname[i]), "CE1") == 0)
143 || (std::strcmp(*(atoms->atomname[i]), "OE1") == 0)
144 || (std::strcmp(*(atoms->atomname[i]), "NE") == 0))
148 else if ((std::strcmp(*(atoms->atomname[i]), "CZ") == 0)
149 || (std::strcmp(*(atoms->atomname[i]), "NZ") == 0))
153 /* HChi flag here too */
154 else if (bHChi && (std::strcmp(*(atoms->atomname[i]), "NH1") == 0))
161 thisres = *(atoms->resinfo[ires].name);
163 /* added by grs - special case for aromatics, whose chis above 2 are
164 not real and produce rubbish output - so set back to -1 */
165 if (std::strcmp(thisres, "PHE") == 0 || std::strcmp(thisres, "TYR") == 0
166 || std::strcmp(thisres, "PTR") == 0 || std::strcmp(thisres, "TRP") == 0
167 || std::strcmp(thisres, "HIS") == 0 || std::strcmp(thisres, "HISA") == 0
168 || std::strcmp(thisres, "HISB") == 0)
170 for (ii = 5; ii <= 7; ii++)
175 /* end fixing aromatics */
177 /* Special case for Pro, has no H */
178 if (std::strcmp(thisres, "PRO") == 0)
182 /* Carbon from previous residue */
187 /* Alpha-carbon from previous residue */
188 if (prev.Cn[1] != -1)
190 atm.minCalpha = prev.Cn[1];
194 /* Check how many dihedrals we have */
195 if ((atm.N != -1) && (atm.Cn[1] != -1) && (atm.C != -1) && (atm.O != -1)
196 && ((atm.H != -1) || (atm.minC != -1)))
198 dl[nl].resnr = ires + 1;
200 dl[nl].atm.Cn[0] = atm.N;
201 if ((atm.Cn[3] != -1) && (atm.Cn[2] != -1) && (atm.Cn[1] != -1))
225 if ((atm.minC != -1) && (atm.minCalpha != -1))
230 /* Prevent use of 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 (residueTypeMap.find(thisres) == residueTypeMap.end())
236 "Unknown residue %s when searching for residue type.\n"
237 "Maybe you need to add a custom residue in residuetypes.dat.",
240 dl[nl].residueName = thisres;
242 sprintf(dl[nl].name, "%s%d", thisres, ires + r0);
248 "Could not find N atom but could find other atoms"
249 " in residue %s%d\n",
254 // Leave only the residues that were recognized to contain dihedrals
257 fprintf(stderr, "\n");
259 fprintf(log, "There are %d residues with dihedrals\n", nl);
271 for (i = 0; (i < maxchi); i++)
276 fprintf(log, "There are %d dihedrals\n", j);
277 fprintf(log, "Dihedral: ");
280 fprintf(log, " Phi ");
284 fprintf(log, " Psi ");
288 for (i = 0; (i < maxchi); i++)
290 fprintf(log, "Chi%d ", i + 1);
293 fprintf(log, "\nNumber: ");
296 fprintf(log, "%4d ", nl);
300 fprintf(log, "%4d ", nl);
304 for (i = 0; (i < maxchi); i++)
306 fprintf(log, "%4d ", nc[i]);
314 gmx_bool has_dihedral(int Dih, const t_dlist& dl)
322 b = ((dl.atm.H != -1) && (dl.atm.N != -1) && (dl.atm.Cn[1] != -1) && (dl.atm.C != -1));
325 b = ((dl.atm.N != -1) && (dl.atm.Cn[1] != -1) && (dl.atm.C != -1) && (dl.atm.O != -1));
328 b = ((dl.atm.minCalpha != -1) && (dl.atm.minC != -1) && (dl.atm.N != -1)
329 && (dl.atm.Cn[1] != -1));
338 b = ((dl.atm.Cn[ddd] != -1) && (dl.atm.Cn[ddd + 1] != -1) && (dl.atm.Cn[ddd + 2] != -1)
339 && (dl.atm.Cn[ddd + 3] != -1));
342 pr_dlist(stdout, gmx::constArrayRefFromArray(&dl, 1), 1, 0, TRUE, TRUE, TRUE, TRUE, MAXCHI);
343 gmx_fatal(FARGS, "Non existent dihedral %d in file %s, line %d", Dih, __FILE__, __LINE__);
348 static void pr_one_ro(FILE* fp, const t_dlist& dl, int nDih, real gmx_unused dt)
351 for (k = 0; k < NROT; k++)
353 fprintf(fp, " %6.2f", dl.rot_occ[nDih][k]);
358 static void pr_ntr_s2(FILE* fp, const t_dlist& dl, int nDih, real dt)
360 fprintf(fp, " %6.2f %6.2f\n", (dt == 0) ? 0 : dl.ntr[nDih] / dt, dl.S2[nDih]);
363 void pr_dlist(FILE* fp,
364 gmx::ArrayRef<const t_dlist> dlist,
373 void (*pr_props)(FILE*, const t_dlist&, int, real);
375 /* Analysis of dihedral transitions etc */
377 if (printtype == edPrintST)
379 pr_props = pr_ntr_s2;
380 fprintf(stderr, "Now printing out transitions and OPs...\n");
384 pr_props = pr_one_ro;
385 fprintf(stderr, "Now printing out rotamer occupancies...\n");
386 fprintf(fp, "\nXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXXX\n\n");
389 /* change atom numbers from 0 based to 1 based */
390 for (const auto& dihedral : dlist)
392 fprintf(fp, "Residue %s\n", dihedral.name);
393 if (printtype == edPrintST)
396 " Angle [ AI, AJ, AK, AL] #tr/ns S^2D \n"
397 "--------------------------------------------\n");
402 " Angle [ AI, AJ, AK, AL] rotamers 0 g(-) t g(+)\n"
403 "--------------------------------------------\n");
408 " Phi [%5d,%5d,%5d,%5d]",
409 (dihedral.atm.H == -1) ? 1 + dihedral.atm.minC : 1 + dihedral.atm.H,
411 1 + dihedral.atm.Cn[1],
413 pr_props(fp, dihedral, edPhi, dt);
418 " Psi [%5d,%5d,%5d,%5d]",
420 1 + dihedral.atm.Cn[1],
423 pr_props(fp, dihedral, edPsi, dt);
425 if (bOmega && has_dihedral(edOmega, dihedral))
428 " Omega [%5d,%5d,%5d,%5d]",
429 1 + dihedral.atm.minCalpha,
430 1 + dihedral.atm.minC,
432 1 + dihedral.atm.Cn[1]);
433 pr_props(fp, dihedral, edOmega, dt);
435 for (int Xi = 0; Xi < MAXCHI; Xi++)
437 if (bChi && (Xi < maxchi) && (dihedral.atm.Cn[Xi + 3] != -1))
440 " Chi%d[%5d,%5d,%5d,%5d]",
442 1 + dihedral.atm.Cn[Xi],
443 1 + dihedral.atm.Cn[Xi + 1],
444 1 + dihedral.atm.Cn[Xi + 2],
445 1 + dihedral.atm.Cn[Xi + 3]);
446 pr_props(fp, dihedral, Xi + edChi1, dt); /* Xi+2 was wrong here */