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,2017,2018,2019, 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.
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.
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.
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.
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.
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.
37 /* This file is completely threadsafe - keep it that way! */
46 #include "gromacs/fileio/pdbio.h"
47 #include "gromacs/gmxpreprocess/pdb2top.h"
48 #include "gromacs/gmxpreprocess/toputil.h"
49 #include "gromacs/math/functions.h"
50 #include "gromacs/math/units.h"
51 #include "gromacs/math/vec.h"
52 #include "gromacs/topology/block.h"
53 #include "gromacs/topology/symtab.h"
54 #include "gromacs/utility/arraysize.h"
55 #include "gromacs/utility/cstringutil.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/smalloc.h"
59 static int in_strings(char* key, int nstr, const char** str)
63 for (j = 0; (j < nstr); j++)
65 if (strcmp(str[j], key) == 0)
74 static bool hbond(rvec x[], int i, int j, real distance)
76 real tol = distance * distance;
79 rvec_sub(x[i], x[j], tmp);
81 return (iprod(tmp, tmp) < tol);
84 static void chk_allhb(t_atoms* pdba, rvec x[], t_blocka* hb, const bool donor[], const bool accept[], real dist)
86 int i, j, k, ii, natom;
89 snew(hb->index, natom + 1);
90 snew(hb->a, 6 * natom);
96 for (i = 0; (i < natom); i++)
100 for (j = i + 1; (j < natom); j++)
102 if ((accept[j]) && (hbond(x, i, j, dist)))
110 for (j = i + 1; (j < natom); j++)
112 if ((donor[j]) && (hbond(x, i, j, dist)))
123 static bool chk_hbonds(int i, t_atoms* pdba, rvec x[], const bool ad[], bool hbond[], rvec xh, real angle, real dist)
126 int j, aj, ri, natom;
132 ri = pdba->atom[i].resind;
133 dist2 = gmx::square(dist);
134 for (j = 0; (j < natom); j++)
136 /* Check whether the other atom is a donor/acceptor and not i */
137 if ((ad[j]) && (j != i))
139 /* Check whether the other atom is on the same ring as well */
140 if ((pdba->atom[j].resind != ri)
141 || ((strcmp(*pdba->atomname[j], "ND1") != 0) && (strcmp(*pdba->atomname[j], "NE2") != 0)))
144 d2 = distance2(x[i], x[j]);
145 rvec_sub(x[i], xh, nh);
146 rvec_sub(x[aj], xh, oh);
147 a = RAD2DEG * acos(cos_angle(nh, oh));
148 if ((d2 < dist2) && (a > angle))
159 static void calc_ringh(rvec xattach, rvec xb, rvec xc, rvec xh)
164 /* Add a proton on a ring to atom attach at distance 0.1 nm */
165 rvec_sub(xattach, xb, tab);
166 rvec_sub(xattach, xc, tac);
167 rvec_add(tab, tac, xh);
170 rvec_inc(xh, xattach);
173 void set_histp(t_atoms* pdba, rvec* x, t_symtab* symtab, real angle, real dist)
175 static const char* prot_acc[] = { "O", "OD1", "OD2", "OE1", "OE2", "OG", "OG1", "OH", "OW" };
176 #define NPA asize(prot_acc)
177 static const char* prot_don[] = { "N", "NH1", "NH2", "NE", "ND1", "ND2", "NE2",
178 "NZ", "OG", "OG1", "OH", "NE1", "OW" };
179 #define NPD asize(prot_don)
181 bool * donor, *acceptor;
186 int i, j, nd, na, hisind, type = -1;
187 int nd1, ne2, cg, cd2, ce1;
194 while (i < natom && gmx_strcasecmp(*pdba->resinfo[pdba->atom[i].resind].name, "HIS") != 0)
203 /* A histidine residue exists that requires automated assignment, so
204 * doing the analysis of donors and acceptors is worthwhile. */
206 "Analysing hydrogen-bonding network for automated assignment of histidine\n"
210 snew(acceptor, natom);
215 for (j = 0; (j < natom); j++)
217 if (in_strings(*pdba->atomname[j], NPA, prot_acc) != -1)
222 if (in_strings(*pdba->atomname[j], NPD, prot_don) != -1)
228 fprintf(stderr, " %d donors and %d acceptors were found.\n", nd, na);
229 chk_allhb(pdba, x, hb, donor, acceptor, dist);
230 fprintf(stderr, "There are %d hydrogen bonds\n", hb->nra);
232 /* Now do the HIS stuff */
236 if (gmx_strcasecmp(*pdba->resinfo[pdba->atom[i].resind].name, "HIS") != 0)
242 if (pdba->atom[i].resind != hisind)
244 hisind = pdba->atom[i].resind;
246 /* Find the atoms in the ring */
247 nd1 = ne2 = cg = cd2 = ce1 = -1;
248 while (i < natom && pdba->atom[i].resind == hisind)
250 atomnm = *pdba->atomname[i];
251 if (strcmp(atomnm, "CD2") == 0)
255 else if (strcmp(atomnm, "CG") == 0)
259 else if (strcmp(atomnm, "CE1") == 0)
263 else if (strcmp(atomnm, "ND1") == 0)
267 else if (strcmp(atomnm, "NE2") == 0)
275 if (!((cg == -1) || (cd2 == -1) || (ce1 == -1) || (nd1 == -1) || (ne2 == -1)))
277 calc_ringh(x[nd1], x[cg], x[ce1], xh1);
278 calc_ringh(x[ne2], x[ce1], x[cd2], xh2);
280 bHDd = chk_hbonds(nd1, pdba, x, acceptor, hbond, xh1, angle, dist);
281 chk_hbonds(nd1, pdba, x, donor, hbond, xh1, angle, dist);
282 bHEd = chk_hbonds(ne2, pdba, x, acceptor, hbond, xh2, angle, dist);
283 chk_hbonds(ne2, pdba, x, donor, hbond, xh2, angle, dist);
300 fprintf(stderr, "Will use %s for residue %d\n", hh[type], pdba->resinfo[hisind].nr);
304 gmx_fatal(FARGS, "Incomplete ring in HIS%d", pdba->resinfo[hisind].nr);
307 pdba->resinfo[hisind].rtp = put_symtab(symtab, hh[type]);