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