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