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