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