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