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 by the GROMACS development team.
7 * Copyright (c) 2019,2020, 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.
38 /* This file is completely threadsafe - keep it that way! */
47 #include "gromacs/fileio/pdbio.h"
48 #include "gromacs/gmxpreprocess/pdb2top.h"
49 #include "gromacs/gmxpreprocess/toputil.h"
50 #include "gromacs/math/functions.h"
51 #include "gromacs/math/units.h"
52 #include "gromacs/math/vec.h"
53 #include "gromacs/topology/block.h"
54 #include "gromacs/topology/symtab.h"
55 #include "gromacs/utility/arraysize.h"
56 #include "gromacs/utility/cstringutil.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/smalloc.h"
60 static int in_strings(char* key
, int nstr
, const char** str
)
64 for (j
= 0; (j
< nstr
); j
++)
66 if (strcmp(str
[j
], key
) == 0)
75 static bool hbond(rvec x
[], int i
, int j
, real distance
)
77 real tol
= distance
* distance
;
80 rvec_sub(x
[i
], x
[j
], tmp
);
82 return (iprod(tmp
, tmp
) < tol
);
85 static void chk_allhb(t_atoms
* pdba
, rvec x
[], t_blocka
* hb
, const bool donor
[], const bool accept
[], real dist
)
87 int i
, j
, k
, ii
, natom
;
90 snew(hb
->index
, natom
+ 1);
91 snew(hb
->a
, 6 * natom
);
97 for (i
= 0; (i
< natom
); i
++)
101 for (j
= i
+ 1; (j
< natom
); j
++)
103 if ((accept
[j
]) && (hbond(x
, i
, j
, dist
)))
111 for (j
= i
+ 1; (j
< natom
); j
++)
113 if ((donor
[j
]) && (hbond(x
, i
, j
, dist
)))
124 static bool chk_hbonds(int i
, t_atoms
* pdba
, rvec x
[], const bool ad
[], bool hbond
[], rvec xh
, real angle
, real dist
)
127 int j
, aj
, ri
, natom
;
133 ri
= pdba
->atom
[i
].resind
;
134 dist2
= gmx::square(dist
);
135 for (j
= 0; (j
< natom
); j
++)
137 /* Check whether the other atom is a donor/acceptor and not i */
138 if ((ad
[j
]) && (j
!= i
))
140 /* Check whether the other atom is on the same ring as well */
141 if ((pdba
->atom
[j
].resind
!= ri
)
142 || ((strcmp(*pdba
->atomname
[j
], "ND1") != 0) && (strcmp(*pdba
->atomname
[j
], "NE2") != 0)))
145 d2
= distance2(x
[i
], x
[j
]);
146 rvec_sub(x
[i
], xh
, nh
);
147 rvec_sub(x
[aj
], xh
, oh
);
148 a
= RAD2DEG
* acos(cos_angle(nh
, oh
));
149 if ((d2
< dist2
) && (a
> angle
))
160 static void calc_ringh(rvec xattach
, rvec xb
, rvec xc
, rvec xh
)
165 /* Add a proton on a ring to atom attach at distance 0.1 nm */
166 rvec_sub(xattach
, xb
, tab
);
167 rvec_sub(xattach
, xc
, tac
);
168 rvec_add(tab
, tac
, xh
);
171 rvec_inc(xh
, xattach
);
174 void set_histp(t_atoms
* pdba
, rvec
* x
, t_symtab
* symtab
, real angle
, real dist
)
176 static const char* prot_acc
[] = { "O", "OD1", "OD2", "OE1", "OE2", "OG", "OG1", "OH", "OW" };
177 #define NPA asize(prot_acc)
178 static const char* prot_don
[] = { "N", "NH1", "NH2", "NE", "ND1", "ND2", "NE2",
179 "NZ", "OG", "OG1", "OH", "NE1", "OW" };
180 #define NPD asize(prot_don)
182 bool * donor
, *acceptor
;
187 int i
, j
, nd
, na
, hisind
, type
= -1;
188 int nd1
, ne2
, cg
, cd2
, ce1
;
195 while (i
< natom
&& gmx_strcasecmp(*pdba
->resinfo
[pdba
->atom
[i
].resind
].name
, "HIS") != 0)
204 /* A histidine residue exists that requires automated assignment, so
205 * doing the analysis of donors and acceptors is worthwhile. */
207 "Analysing hydrogen-bonding network for automated assignment of histidine\n"
211 snew(acceptor
, natom
);
216 for (j
= 0; (j
< natom
); j
++)
218 if (in_strings(*pdba
->atomname
[j
], NPA
, prot_acc
) != -1)
223 if (in_strings(*pdba
->atomname
[j
], NPD
, prot_don
) != -1)
229 fprintf(stderr
, " %d donors and %d acceptors were found.\n", nd
, na
);
230 chk_allhb(pdba
, x
, hb
, donor
, acceptor
, dist
);
231 fprintf(stderr
, "There are %d hydrogen bonds\n", hb
->nra
);
233 /* Now do the HIS stuff */
237 if (gmx_strcasecmp(*pdba
->resinfo
[pdba
->atom
[i
].resind
].name
, "HIS") != 0)
243 if (pdba
->atom
[i
].resind
!= hisind
)
245 hisind
= pdba
->atom
[i
].resind
;
247 /* Find the atoms in the ring */
248 nd1
= ne2
= cg
= cd2
= ce1
= -1;
249 while (i
< natom
&& pdba
->atom
[i
].resind
== hisind
)
251 atomnm
= *pdba
->atomname
[i
];
252 if (strcmp(atomnm
, "CD2") == 0)
256 else if (strcmp(atomnm
, "CG") == 0)
260 else if (strcmp(atomnm
, "CE1") == 0)
264 else if (strcmp(atomnm
, "ND1") == 0)
268 else if (strcmp(atomnm
, "NE2") == 0)
276 if (!((cg
== -1) || (cd2
== -1) || (ce1
== -1) || (nd1
== -1) || (ne2
== -1)))
278 calc_ringh(x
[nd1
], x
[cg
], x
[ce1
], xh1
);
279 calc_ringh(x
[ne2
], x
[ce1
], x
[cd2
], xh2
);
281 bHDd
= chk_hbonds(nd1
, pdba
, x
, acceptor
, hbond
, xh1
, angle
, dist
);
282 chk_hbonds(nd1
, pdba
, x
, donor
, hbond
, xh1
, angle
, dist
);
283 bHEd
= chk_hbonds(ne2
, pdba
, x
, acceptor
, hbond
, xh2
, angle
, dist
);
284 chk_hbonds(ne2
, pdba
, x
, donor
, hbond
, xh2
, angle
, dist
);
301 fprintf(stderr
, "Will use %s for residue %d\n", hh
[type
], pdba
->resinfo
[hisind
].nr
);
305 gmx_fatal(FARGS
, "Incomplete ring in HIS%d", pdba
->resinfo
[hisind
].nr
);
308 pdba
->resinfo
[hisind
].rtp
= put_symtab(symtab
, hh
[type
]);