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