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/math/vec.h"
47 #include "gromacs/math/units.h"
50 #include "gromacs/utility/cstringutil.h"
53 #include "gromacs/topology/block.h"
54 #include "gromacs/utility/fatalerror.h"
55 #include "gromacs/utility/smalloc.h"
57 static int in_strings(char *key, int nstr, const char **str)
61 for (j = 0; (j < nstr); j++)
63 if (strcmp(str[j], key) == 0)
72 static gmx_bool hbond(rvec x[], int i, int j, real distance)
74 real tol = distance*distance;
77 rvec_sub(x[i], x[j], tmp);
79 return (iprod(tmp, tmp) < tol);
82 static void chk_allhb(t_atoms *pdba, rvec x[], t_blocka *hb,
83 gmx_bool donor[], gmx_bool accept[], real dist)
85 int i, j, k, ii, natom;
88 snew(hb->index, natom+1);
95 for (i = 0; (i < natom); i++)
99 for (j = i+1; (j < natom); j++)
101 if ((accept[j]) && (hbond(x, i, j, dist)))
109 for (j = i+1; (j < natom); j++)
111 if ((donor[j]) && (hbond(x, i, j, dist)))
122 static void pr_hbonds(FILE *fp, t_blocka *hb, t_atoms *pdba)
126 fprintf(fp, "Dumping all hydrogen bonds!\n");
127 for (i = 0; (i < hb->nr); i++)
131 for (j = j0; (j < j1); j++)
134 fprintf(fp, "%5s%4d%5s - %5s%4d%5s\n",
135 *pdba->resinfo[pdba->atom[i].resind].name,
136 pdba->resinfo[pdba->atom[i].resind].nr, *pdba->atomname[i],
137 *pdba->resinfo[pdba->atom[k].resind].name,
138 pdba->resinfo[pdba->atom[k].resind].nr, *pdba->atomname[k]);
143 static gmx_bool chk_hbonds(int i, t_atoms *pdba, rvec x[],
144 gmx_bool ad[], gmx_bool hbond[], rvec xh,
145 real angle, real dist)
148 int j, aj, ri, natom;
154 ri = pdba->atom[i].resind;
156 for (j = 0; (j < natom); j++)
158 /* Check whether the other atom is a donor/acceptor and not i */
159 if ((ad[j]) && (j != i))
161 /* Check whether the other atom is on the same ring as well */
162 if ((pdba->atom[j].resind != ri) ||
163 ((strcmp(*pdba->atomname[j], "ND1") != 0) &&
164 (strcmp(*pdba->atomname[j], "NE2") != 0)))
167 d2 = distance2(x[i], x[j]);
168 rvec_sub(x[i], xh, nh);
169 rvec_sub(x[aj], xh, oh);
170 a = RAD2DEG * acos(cos_angle(nh, oh));
171 if ((d2 < dist2) && (a > angle))
176 "HBOND between %s%d-%s and %s%d-%s is %g nm, %g deg\n",
177 *pdba->resinfo[pdba->atom[i].resind].name,
178 pdba->resinfo[pdba->atom[i].resind].nr, *pdba->atomname[i],
179 *pdba->resinfo[pdba->atom[aj].resind].name,
180 pdba->resinfo[pdba->atom[aj].resind].nr, *pdba->atomname[aj],
192 static void calc_ringh(rvec xattach, rvec xb, rvec xc, rvec xh)
197 /* Add a proton on a ring to atom attach at distance 0.1 nm */
198 rvec_sub(xattach, xb, tab);
199 rvec_sub(xattach, xc, tac);
200 rvec_add(tab, tac, xh);
203 rvec_inc(xh, xattach);
206 void set_histp(t_atoms *pdba, rvec *x, real angle, real dist)
208 static const char *prot_acc[] = {
209 "O", "OD1", "OD2", "OE1", "OE2", "OG", "OG1", "OH", "OW"
211 #define NPA asize(prot_acc)
212 static const char *prot_don[] = {
213 "N", "NH1", "NH2", "NE", "ND1", "ND2", "NE2", "NZ", "OG", "OG1", "OH", "NE1", "OW"
215 #define NPD asize(prot_don)
217 gmx_bool *donor, *acceptor;
218 gmx_bool *hbond, bHaveH = FALSE;
222 int i, j, nd, na, aj, hisind, his0, type = -1;
223 int nd1, ne2, cg, cd2, ce1;
232 gmx_strcasecmp(*pdba->resinfo[pdba->atom[i].resind].name, "HIS") != 0)
241 /* A histidine residue exists that requires automated assignment, so
242 * doing the analysis of donors and acceptors is worthwhile. */
244 "Analysing hydrogen-bonding network for automated assignment of histidine\n"
248 snew(acceptor, natom);
253 for (j = 0; (j < natom); j++)
255 if (in_strings(*pdba->atomname[j], NPA, prot_acc) != -1)
260 if (in_strings(*pdba->atomname[j], NPD, prot_don) != -1)
266 fprintf(stderr, " %d donors and %d acceptors were found.\n", nd, na);
267 chk_allhb(pdba, x, hb, donor, acceptor, dist);
270 pr_hbonds(debug, hb, pdba);
272 fprintf(stderr, "There are %d hydrogen bonds\n", hb->nra);
274 /* Now do the HIS stuff */
278 if (gmx_strcasecmp(*pdba->resinfo[pdba->atom[i].resind].name, "HIS") != 0)
284 if (pdba->atom[i].resind != hisind)
286 hisind = pdba->atom[i].resind;
288 /* Find the atoms in the ring */
289 nd1 = ne2 = cg = cd2 = ce1 = -1;
290 while (i < natom && pdba->atom[i].resind == hisind)
292 atomnm = *pdba->atomname[i];
293 if (strcmp(atomnm, "CD2") == 0)
297 else if (strcmp(atomnm, "CG") == 0)
301 else if (strcmp(atomnm, "CE1") == 0)
305 else if (strcmp(atomnm, "ND1") == 0)
309 else if (strcmp(atomnm, "NE2") == 0)
317 if (!((cg == -1 ) || (cd2 == -1) || (ce1 == -1) ||
318 (nd1 == -1) || (ne2 == -1)))
320 calc_ringh(x[nd1], x[cg], x[ce1], xh1);
321 calc_ringh(x[ne2], x[ce1], x[cd2], xh2);
323 bHDd = chk_hbonds(nd1, pdba, x, acceptor, hbond, xh1, angle, dist);
324 chk_hbonds(nd1, pdba, x, donor, hbond, xh1, angle, dist);
325 bHEd = chk_hbonds(ne2, pdba, x, acceptor, hbond, xh2, angle, dist);
326 chk_hbonds(ne2, pdba, x, donor, hbond, xh2, angle, dist);
343 fprintf(stderr, "Will use %s for residue %d\n",
344 hh[type], pdba->resinfo[hisind].nr);
348 gmx_fatal(FARGS, "Incomplete ring in HIS%d",
349 pdba->resinfo[hisind].nr);
352 snew(pdba->resinfo[hisind].rtp, 1);
353 *pdba->resinfo[hisind].rtp = strdup(hh[type]);