Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / specbond.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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,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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #include "gmxpre.h"
38
39 #include "specbond.h"
40
41 #include <cctype>
42 #include <cmath>
43 #include <cstring>
44
45 #include <algorithm>
46
47 #include "gromacs/fileio/pdbio.h"
48 #include "gromacs/gmxpreprocess/pdb2top.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/utility/cstringutil.h"
51 #include "gromacs/utility/fatalerror.h"
52 #include "gromacs/utility/smalloc.h"
53 #include "gromacs/utility/strdb.h"
54
55 struct SpecialBond
56 {
57     std::string firstResidue, secondResidue;
58     std::string firstAtomName, secondAtomName;
59     std::string newFirstResidue, newSecondResidue;
60     int         firstBondNumber, secondBondNumber;
61     real        length;
62 };
63
64 static bool yesno()
65 {
66     char c;
67
68     do
69     {
70         c = toupper(fgetc(stdin));
71     } while ((c != 'Y') && (c != 'N'));
72
73     return (c == 'Y');
74 }
75
76 std::vector<SpecialBond> generateSpecialBonds()
77 {
78     const char* sbfile = "specbond.dat";
79
80     std::vector<SpecialBond> specialBonds;
81     char                     r1buf[32], r2buf[32], a1buf[32], a2buf[32], nr1buf[32], nr2buf[32];
82     double                   length;
83     int                      nb1, nb2;
84     char**                   lines;
85
86     int nlines = get_lines(sbfile, &lines);
87     for (int i = 0; (i < nlines); i++)
88     {
89         if (sscanf(lines[i], "%s%s%d%s%s%d%lf%s%s", r1buf, a1buf, &nb1, r2buf, a2buf, &nb2, &length,
90                    nr1buf, nr2buf)
91             != 9)
92         {
93             fprintf(stderr, "Invalid line '%s' in %s\n", lines[i], sbfile);
94         }
95         else
96         {
97             SpecialBond newBond;
98
99             newBond.firstResidue     = r1buf;
100             newBond.secondResidue    = r2buf;
101             newBond.newFirstResidue  = nr1buf;
102             newBond.newSecondResidue = nr2buf;
103             newBond.firstAtomName    = a1buf;
104             newBond.secondAtomName   = a2buf;
105             newBond.firstBondNumber  = nb1;
106             newBond.secondBondNumber = nb2;
107             newBond.length           = length;
108             specialBonds.push_back(newBond);
109         }
110         sfree(lines[i]);
111     }
112     if (nlines > 0)
113     {
114         sfree(lines);
115     }
116     fprintf(stderr, "%zu out of %d lines of %s converted successfully\n", specialBonds.size(),
117             nlines, sbfile);
118
119     return specialBonds;
120 }
121
122 static bool is_special(gmx::ArrayRef<const SpecialBond> sb, const char* res, const char* atom)
123 {
124     for (const auto& bond : sb)
125     {
126         if (((strncmp(bond.firstResidue.c_str(), res, 3) == 0)
127              && (gmx::equalCaseInsensitive(bond.firstAtomName, atom)))
128             || ((strncmp(bond.secondResidue.c_str(), res, 3) == 0)
129                 && (gmx::equalCaseInsensitive(bond.secondAtomName, atom))))
130         {
131             return TRUE;
132         }
133     }
134     return FALSE;
135 }
136
137 static bool is_bond(gmx::ArrayRef<const SpecialBond> sb, t_atoms* pdba, int a1, int a2, real d, int* index_sb, bool* bSwap)
138 {
139     const char* at1  = *pdba->atomname[a1];
140     const char* at2  = *pdba->atomname[a2];
141     const char* res1 = *pdba->resinfo[pdba->atom[a1].resind].name;
142     const char* res2 = *pdba->resinfo[pdba->atom[a2].resind].name;
143
144     int i = 0;
145     for (const auto& bond : sb)
146     {
147         *index_sb = i;
148         if (((strncmp(bond.firstResidue.c_str(), res1, 3) == 0)
149              && (gmx::equalCaseInsensitive(bond.firstAtomName, at1))
150              && (strncmp(bond.secondResidue.c_str(), res2, 3) == 0)
151              && (gmx::equalCaseInsensitive(bond.secondAtomName, at2))))
152         {
153             *bSwap = FALSE;
154             if ((0.9 * bond.length < d) && (1.1 * bond.length > d))
155             {
156                 return TRUE;
157             }
158         }
159         if (((strncmp(bond.firstResidue.c_str(), res2, 3) == 0)
160              && (gmx::equalCaseInsensitive(bond.firstAtomName, at2))
161              && (strncmp(bond.secondResidue.c_str(), res1, 3) == 0)
162              && (gmx::equalCaseInsensitive(bond.secondAtomName, at1))))
163         {
164             *bSwap = TRUE;
165             if ((0.9 * bond.length < d) && (1.1 * bond.length > d))
166             {
167                 return TRUE;
168             }
169         }
170         i++;
171     }
172     return FALSE;
173 }
174
175 static void rename_1res(t_atoms* pdba, int resind, const char* newres, bool bVerbose)
176 {
177     if (bVerbose)
178     {
179         printf("Using rtp entry %s for %s %d\n", newres, *pdba->resinfo[resind].name,
180                pdba->resinfo[resind].nr);
181     }
182     /* this used to free *resname, which messes up the symtab! */
183     snew(pdba->resinfo[resind].rtp, 1);
184     *pdba->resinfo[resind].rtp = gmx_strdup(newres);
185 }
186
187 std::vector<DisulfideBond> makeDisulfideBonds(t_atoms* pdba, rvec x[], bool bInteractive, bool bVerbose)
188 {
189     std::vector<SpecialBond>   specialBonds = generateSpecialBonds();
190     std::vector<DisulfideBond> bonds;
191     bool                       bSwap;
192     int                        index_sb;
193     char                       buf[10];
194
195
196     if (!specialBonds.empty())
197     {
198         std::vector<int> specialBondResIdxs;
199         std::vector<int> specialBondAtomIdxs;
200
201         for (int i = 0; (i < pdba->nr); i++)
202         {
203             /* Check if this atom is special and if it is not a double atom
204              * in the input that still needs to be removed.
205              */
206             int prevAtom = -1;
207             if (!specialBondAtomIdxs.empty())
208             {
209                 prevAtom = specialBondAtomIdxs.back();
210             }
211             if (is_special(specialBonds, *pdba->resinfo[pdba->atom[i].resind].name, *pdba->atomname[i])
212                 && !(!specialBondAtomIdxs.empty() && pdba->atom[prevAtom].resind == pdba->atom[i].resind
213                      && gmx_strcasecmp(*pdba->atomname[prevAtom], *pdba->atomname[i]) == 0))
214             {
215                 specialBondResIdxs.push_back(pdba->atom[i].resind);
216                 specialBondAtomIdxs.push_back(i);
217             }
218         }
219         /* distance matrix d[nspec][nspec] */
220         int                            nspec = specialBondAtomIdxs.size();
221         std::vector<std::vector<real>> d(nspec);
222         for (int i = 0; (i < nspec); i++)
223         {
224             d[i].resize(nspec);
225             int ai = specialBondAtomIdxs[i];
226             for (int j = 0; (j < nspec); j++)
227             {
228                 int aj  = specialBondAtomIdxs[j];
229                 d[i][j] = std::sqrt(distance2(x[ai], x[aj]));
230             }
231         }
232         if (nspec > 1)
233         {
234 #define MAXCOL 7
235             fprintf(stderr, "Special Atom Distance matrix:\n");
236             for (int b = 0; (b < nspec); b += MAXCOL)
237             {
238                 /* print resname/number column headings */
239                 fprintf(stderr, "%8s%8s", "", "");
240                 int e = std::min(b + MAXCOL, nspec - 1);
241                 for (int i = b; (i < e); i++)
242                 {
243                     sprintf(buf, "%s%d", *pdba->resinfo[pdba->atom[specialBondAtomIdxs[i]].resind].name,
244                             pdba->resinfo[specialBondResIdxs[i]].nr);
245                     fprintf(stderr, "%8s", buf);
246                 }
247                 fprintf(stderr, "\n");
248                 /* print atomname/number column headings */
249                 fprintf(stderr, "%8s%8s", "", "");
250                 e = std::min(b + MAXCOL, nspec - 1);
251                 for (int i = b; (i < e); i++)
252                 {
253                     std::string buf = gmx::formatString("%s%d", *pdba->atomname[specialBondAtomIdxs[i]],
254                                                         specialBondAtomIdxs[i] + 1);
255                     fprintf(stderr, "%8s", buf.c_str());
256                 }
257                 fprintf(stderr, "\n");
258                 /* print matrix */
259                 e = std::min(b + MAXCOL, nspec);
260                 for (int i = b + 1; (i < nspec); i++)
261                 {
262                     std::string buf = gmx::formatString(
263                             "%s%d", *pdba->resinfo[pdba->atom[specialBondAtomIdxs[i]].resind].name,
264                             pdba->resinfo[specialBondResIdxs[i]].nr);
265                     fprintf(stderr, "%8s", buf.c_str());
266                     buf = gmx::formatString("%s%d", *pdba->atomname[specialBondAtomIdxs[i]],
267                                             specialBondAtomIdxs[i] + 1);
268                     fprintf(stderr, "%8s", buf.c_str());
269                     int e2 = std::min(i, e);
270                     for (int j = b; (j < e2); j++)
271                     {
272                         fprintf(stderr, " %7.3f", d[i][j]);
273                     }
274                     fprintf(stderr, "\n");
275                 }
276             }
277         }
278
279         for (int i = 0; (i < nspec); i++)
280         {
281             int ai = specialBondAtomIdxs[i];
282             for (int j = i + 1; (j < nspec); j++)
283             {
284                 int aj = specialBondAtomIdxs[j];
285                 /* Ensure creation of at most nspec special bonds to avoid overflowing bonds[] */
286                 if (bonds.size() < specialBondAtomIdxs.size()
287                     && is_bond(specialBonds, pdba, ai, aj, d[i][j], &index_sb, &bSwap))
288                 {
289                     fprintf(stderr, "%s %s-%d %s-%d and %s-%d %s-%d%s", bInteractive ? "Link" : "Linking",
290                             *pdba->resinfo[pdba->atom[ai].resind].name,
291                             pdba->resinfo[specialBondResIdxs[i]].nr, *pdba->atomname[ai], ai + 1,
292                             *pdba->resinfo[pdba->atom[aj].resind].name,
293                             pdba->resinfo[specialBondResIdxs[j]].nr, *pdba->atomname[aj], aj + 1,
294                             bInteractive ? " (y/n) ?" : "...\n");
295                     bool bDoit = bInteractive ? yesno() : true;
296
297                     if (bDoit)
298                     {
299                         DisulfideBond newBond;
300                         /* Store the residue numbers in the bonds array */
301                         newBond.firstResidue  = specialBondResIdxs[i];
302                         newBond.secondResidue = specialBondResIdxs[j];
303                         newBond.firstAtom     = *pdba->atomname[ai];
304                         newBond.secondAtom    = *pdba->atomname[aj];
305                         bonds.push_back(newBond);
306                         /* rename residues */
307                         if (bSwap)
308                         {
309                             rename_1res(pdba, specialBondResIdxs[i],
310                                         specialBonds[index_sb].newSecondResidue.c_str(), bVerbose);
311                             rename_1res(pdba, specialBondResIdxs[j],
312                                         specialBonds[index_sb].newFirstResidue.c_str(), bVerbose);
313                         }
314                         else
315                         {
316                             rename_1res(pdba, specialBondResIdxs[i],
317                                         specialBonds[index_sb].newFirstResidue.c_str(), bVerbose);
318                             rename_1res(pdba, specialBondResIdxs[j],
319                                         specialBonds[index_sb].newSecondResidue.c_str(), bVerbose);
320                         }
321                     }
322                 }
323             }
324         }
325     }
326
327     return bonds;
328 }