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, 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! */
45 #include "gromacs/fileio/pdbio.h"
46 #include "gromacs/utility/smalloc.h"
51 #include "gromacs/utility/cstringutil.h"
54 #include "gromacs/utility/fatalerror.h"
56 static int in_strings(char *key, int nstr, const char **str)
60 for (j = 0; (j < nstr); j++)
62 if (strcmp(str[j], key) == 0)
71 static gmx_bool hbond(rvec x[], int i, int j, real distance)
73 real tol = distance*distance;
76 rvec_sub(x[i], x[j], tmp);
78 return (iprod(tmp, tmp) < tol);
81 static void chk_allhb(t_atoms *pdba, rvec x[], t_blocka *hb,
82 gmx_bool donor[], gmx_bool accept[], real dist)
84 int i, j, k, ii, natom;
87 snew(hb->index, natom+1);
94 for (i = 0; (i < natom); i++)
98 for (j = i+1; (j < natom); j++)
100 if ((accept[j]) && (hbond(x, i, j, dist)))
108 for (j = i+1; (j < natom); j++)
110 if ((donor[j]) && (hbond(x, i, j, dist)))
121 static void pr_hbonds(FILE *fp, t_blocka *hb, t_atoms *pdba)
125 fprintf(fp, "Dumping all hydrogen bonds!\n");
126 for (i = 0; (i < hb->nr); i++)
130 for (j = j0; (j < j1); j++)
133 fprintf(fp, "%5s%4d%5s - %5s%4d%5s\n",
134 *pdba->resinfo[pdba->atom[i].resind].name,
135 pdba->resinfo[pdba->atom[i].resind].nr, *pdba->atomname[i],
136 *pdba->resinfo[pdba->atom[k].resind].name,
137 pdba->resinfo[pdba->atom[k].resind].nr, *pdba->atomname[k]);
142 static gmx_bool chk_hbonds(int i, t_atoms *pdba, rvec x[],
143 gmx_bool ad[], gmx_bool hbond[], rvec xh,
144 real angle, real dist)
147 int j, aj, ri, natom;
153 ri = pdba->atom[i].resind;
155 for (j = 0; (j < natom); j++)
157 /* Check whether the other atom is a donor/acceptor and not i */
158 if ((ad[j]) && (j != i))
160 /* Check whether the other atom is on the same ring as well */
161 if ((pdba->atom[j].resind != ri) ||
162 ((strcmp(*pdba->atomname[j], "ND1") != 0) &&
163 (strcmp(*pdba->atomname[j], "NE2") != 0)))
166 d2 = distance2(x[i], x[j]);
167 rvec_sub(x[i], xh, nh);
168 rvec_sub(x[aj], xh, oh);
169 a = RAD2DEG * acos(cos_angle(nh, oh));
170 if ((d2 < dist2) && (a > angle))
175 "HBOND between %s%d-%s and %s%d-%s is %g nm, %g deg\n",
176 *pdba->resinfo[pdba->atom[i].resind].name,
177 pdba->resinfo[pdba->atom[i].resind].nr, *pdba->atomname[i],
178 *pdba->resinfo[pdba->atom[aj].resind].name,
179 pdba->resinfo[pdba->atom[aj].resind].nr, *pdba->atomname[aj],
191 static void calc_ringh(rvec xattach, rvec xb, rvec xc, rvec xh)
196 /* Add a proton on a ring to atom attach at distance 0.1 nm */
197 rvec_sub(xattach, xb, tab);
198 rvec_sub(xattach, xc, tac);
199 rvec_add(tab, tac, xh);
202 rvec_inc(xh, xattach);
205 void set_histp(t_atoms *pdba, rvec *x, real angle, real dist)
207 static const char *prot_acc[] = {
208 "O", "OD1", "OD2", "OE1", "OE2", "OG", "OG1", "OH", "OW"
210 #define NPA asize(prot_acc)
211 static const char *prot_don[] = {
212 "N", "NH1", "NH2", "NE", "ND1", "ND2", "NE2", "NZ", "OG", "OG1", "OH", "NE1", "OW"
214 #define NPD asize(prot_don)
216 gmx_bool *donor, *acceptor;
217 gmx_bool *hbond, bHaveH = FALSE;
221 int i, j, nd, na, aj, hisind, his0, type = -1;
222 int nd1, ne2, cg, cd2, ce1;
231 gmx_strcasecmp(*pdba->resinfo[pdba->atom[i].resind].name, "HIS") != 0)
240 /* A histidine residue exists that requires automated assignment, so
241 * doing the analysis of donors and acceptors is worthwhile. */
243 "Analysing hydrogen-bonding network for automated assignment of histidine\n"
247 snew(acceptor, natom);
252 for (j = 0; (j < natom); j++)
254 if (in_strings(*pdba->atomname[j], NPA, prot_acc) != -1)
259 if (in_strings(*pdba->atomname[j], NPD, prot_don) != -1)
265 fprintf(stderr, " %d donors and %d acceptors were found.\n", nd, na);
266 chk_allhb(pdba, x, hb, donor, acceptor, dist);
269 pr_hbonds(debug, hb, pdba);
271 fprintf(stderr, "There are %d hydrogen bonds\n", hb->nra);
273 /* Now do the HIS stuff */
277 if (gmx_strcasecmp(*pdba->resinfo[pdba->atom[i].resind].name, "HIS") != 0)
283 if (pdba->atom[i].resind != hisind)
285 hisind = pdba->atom[i].resind;
287 /* Find the atoms in the ring */
288 nd1 = ne2 = cg = cd2 = ce1 = -1;
289 while (i < natom && pdba->atom[i].resind == hisind)
291 atomnm = *pdba->atomname[i];
292 if (strcmp(atomnm, "CD2") == 0)
296 else if (strcmp(atomnm, "CG") == 0)
300 else if (strcmp(atomnm, "CE1") == 0)
304 else if (strcmp(atomnm, "ND1") == 0)
308 else if (strcmp(atomnm, "NE2") == 0)
316 if (!((cg == -1 ) || (cd2 == -1) || (ce1 == -1) ||
317 (nd1 == -1) || (ne2 == -1)))
319 calc_ringh(x[nd1], x[cg], x[ce1], xh1);
320 calc_ringh(x[ne2], x[ce1], x[cd2], xh2);
322 bHDd = chk_hbonds(nd1, pdba, x, acceptor, hbond, xh1, angle, dist);
323 chk_hbonds(nd1, pdba, x, donor, hbond, xh1, angle, dist);
324 bHEd = chk_hbonds(ne2, pdba, x, acceptor, hbond, xh2, angle, dist);
325 chk_hbonds(ne2, pdba, x, donor, hbond, xh2, angle, dist);
342 fprintf(stderr, "Will use %s for residue %d\n",
343 hh[type], pdba->resinfo[hisind].nr);
347 gmx_fatal(FARGS, "Incomplete ring in HIS%d",
348 pdba->resinfo[hisind].nr);
351 snew(pdba->resinfo[hisind].rtp, 1);
352 *pdba->resinfo[hisind].rtp = strdup(hh[type]);