Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / xlate.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,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 "xlate.h"
40
41 #include <cctype>
42 #include <cstring>
43
44 #include <string>
45 #include <vector>
46
47 #include "gromacs/gmxpreprocess/fflibutil.h"
48 #include "gromacs/gmxpreprocess/grompp_impl.h"
49 #include "gromacs/topology/residuetypes.h"
50 #include "gromacs/topology/symtab.h"
51 #include "gromacs/utility/cstringutil.h"
52 #include "gromacs/utility/fatalerror.h"
53 #include "gromacs/utility/futil.h"
54 #include "gromacs/utility/smalloc.h"
55 #include "gromacs/utility/strdb.h"
56
57 #include "hackblock.h"
58
59 typedef struct
60 {
61     char* filebase;
62     char* res;
63     char* atom;
64     char* replace;
65 } t_xlate_atom;
66
67 static void get_xlatoms(const std::string& filename, FILE* fp, int* nptr, t_xlate_atom** xlptr)
68 {
69     char          filebase[STRLEN];
70     char          line[STRLEN];
71     char          abuf[1024], rbuf[1024], repbuf[1024], dumbuf[1024];
72     char*         _ptr;
73     int           n, na, idum;
74     t_xlate_atom* xl;
75
76     fflib_filename_base(filename.c_str(), filebase, STRLEN);
77
78     n  = *nptr;
79     xl = *xlptr;
80
81     while (get_a_line(fp, line, STRLEN))
82     {
83         na = sscanf(line, "%s%s%s%s", rbuf, abuf, repbuf, dumbuf);
84         /* Check if we are reading an old format file with the number of items
85          * on the first line.
86          */
87         if (na == 1 && n == *nptr && sscanf(rbuf, "%d", &idum) == 1)
88         {
89             continue;
90         }
91         if (na != 3)
92         {
93             gmx_fatal(FARGS, "Expected a residue name and two atom names in file '%s', not '%s'",
94                       filename.c_str(), line);
95         }
96
97         srenew(xl, n + 1);
98         xl[n].filebase = gmx_strdup(filebase);
99
100         /* Use wildcards... */
101         if (strcmp(rbuf, "*") != 0)
102         {
103             xl[n].res = gmx_strdup(rbuf);
104         }
105         else
106         {
107             xl[n].res = nullptr;
108         }
109
110         /* Replace underscores in the string by spaces */
111         while ((_ptr = strchr(abuf, '_')) != nullptr)
112         {
113             *_ptr = ' ';
114         }
115
116         xl[n].atom    = gmx_strdup(abuf);
117         xl[n].replace = gmx_strdup(repbuf);
118         n++;
119     }
120
121     *nptr  = n;
122     *xlptr = xl;
123 }
124
125 static void done_xlatom(int nxlate, t_xlate_atom* xlatom)
126 {
127     int i;
128
129     for (i = 0; (i < nxlate); i++)
130     {
131         sfree(xlatom[i].filebase);
132         if (xlatom[i].res != nullptr)
133         {
134             sfree(xlatom[i].res);
135         }
136         sfree(xlatom[i].atom);
137         sfree(xlatom[i].replace);
138     }
139     sfree(xlatom);
140 }
141
142 void rename_atoms(const char*                            xlfile,
143                   const char*                            ffdir,
144                   t_atoms*                               atoms,
145                   t_symtab*                              symtab,
146                   gmx::ArrayRef<const PreprocessResidue> localPpResidue,
147                   bool                                   bResname,
148                   ResidueType*                           rt,
149                   bool                                   bReorderNum,
150                   bool                                   bVerbose)
151 {
152     int           nxlate, a, i, resind;
153     t_xlate_atom* xlatom;
154     char          c, *rnm, atombuf[32];
155     bool          bReorderedNum, bRenamed, bMatch;
156     bool          bStartTerm, bEndTerm;
157
158     nxlate = 0;
159     xlatom = nullptr;
160     if (xlfile != nullptr)
161     {
162         gmx::FilePtr fp = gmx::openLibraryFile(xlfile);
163         get_xlatoms(xlfile, fp.get(), &nxlate, &xlatom);
164     }
165     else
166     {
167         std::vector<std::string> fns = fflib_search_file_end(ffdir, ".arn", FALSE);
168         for (const auto& filename : fns)
169         {
170             FILE* fp = fflib_open(filename);
171             get_xlatoms(filename, fp, &nxlate, &xlatom);
172             gmx_ffclose(fp);
173         }
174     }
175
176     for (a = 0; (a < atoms->nr); a++)
177     {
178         resind = atoms->atom[a].resind;
179
180         bStartTerm = (resind == 0) || atoms->resinfo[resind].chainnum != atoms->resinfo[resind - 1].chainnum;
181         bEndTerm = (resind >= atoms->nres - 1)
182                    || atoms->resinfo[resind].chainnum != atoms->resinfo[resind + 1].chainnum;
183
184         if (bResname)
185         {
186             rnm = *(atoms->resinfo[resind].name);
187         }
188         else
189         {
190             rnm = *(atoms->resinfo[resind].rtp);
191         }
192
193         strcpy(atombuf, *(atoms->atomname[a]));
194         bReorderedNum = FALSE;
195         if (bReorderNum)
196         {
197             if (isdigit(atombuf[0]))
198             {
199                 c = atombuf[0];
200                 for (i = 0; (static_cast<size_t>(i) < strlen(atombuf) - 1); i++)
201                 {
202                     atombuf[i] = atombuf[i + 1];
203                 }
204                 atombuf[i]    = c;
205                 bReorderedNum = TRUE;
206             }
207         }
208         bRenamed = FALSE;
209         for (i = 0; (i < nxlate) && !bRenamed; i++)
210         {
211             /* Check if the base file name of the rtp and arn entry match */
212             if (localPpResidue.empty()
213                 || gmx::equalCaseInsensitive(localPpResidue[resind].filebase, xlatom[i].filebase))
214             {
215                 /* Match the residue name */
216                 bMatch = (xlatom[i].res == nullptr
217                           || (gmx_strcasecmp("protein-nterm", xlatom[i].res) == 0
218                               && rt->namedResidueHasType(rnm, "Protein") && bStartTerm)
219                           || (gmx_strcasecmp("protein-cterm", xlatom[i].res) == 0
220                               && rt->namedResidueHasType(rnm, "Protein") && bEndTerm)
221                           || (gmx_strcasecmp("protein", xlatom[i].res) == 0
222                               && rt->namedResidueHasType(rnm, "Protein"))
223                           || (gmx_strcasecmp("DNA", xlatom[i].res) == 0
224                               && rt->namedResidueHasType(rnm, "DNA"))
225                           || (gmx_strcasecmp("RNA", xlatom[i].res) == 0
226                               && rt->namedResidueHasType(rnm, "RNA")));
227                 if (!bMatch)
228                 {
229                     const char* ptr0 = rnm;
230                     const char* ptr1 = xlatom[i].res;
231                     while (ptr0[0] != '\0' && ptr1[0] != '\0' && (ptr0[0] == ptr1[0] || ptr1[0] == '?'))
232                     {
233                         ptr0++;
234                         ptr1++;
235                     }
236                     bMatch = (ptr0[0] == '\0' && ptr1[0] == '\0');
237                 }
238                 if (bMatch && strcmp(atombuf, xlatom[i].atom) == 0)
239                 {
240                     /* We have a match. */
241                     /* Don't free the old atomname,
242                      * since it might be in the symtab.
243                      */
244                     const char* ptr0 = xlatom[i].replace;
245                     if (bVerbose)
246                     {
247                         printf("Renaming atom '%s' in residue %d %s to '%s'\n", *atoms->atomname[a],
248                                atoms->resinfo[resind].nr, *atoms->resinfo[resind].name, ptr0);
249                     }
250                     atoms->atomname[a] = put_symtab(symtab, ptr0);
251                     bRenamed           = TRUE;
252                 }
253             }
254         }
255         if (bReorderedNum && !bRenamed)
256         {
257             atoms->atomname[a] = put_symtab(symtab, atombuf);
258         }
259     }
260
261     done_xlatom(nxlate, xlatom);
262 }