Adjust more copyright headers
[alexxy/gromacs.git] / src / programs / gmx / 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, 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 "strdb.h"
45 #include "string2.h"
46 #include "smalloc.h"
47 #include "symtab.h"
48 #include "index.h"
49 #include "gromacs/fileio/futil.h"
50 #include "fflibutil.h"
51 #include "hackblock.h"
52 #include "gmx_fatal.h"
53 #include "xlate.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 = strdup(filebase);
94
95         /* Use wildcards... */
96         if (strcmp(rbuf, "*") != 0)
97         {
98             xl[n].res = 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    = strdup(abuf);
112         xl[n].replace = 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
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             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         if (bResname)
175         {
176             rnm = *(atoms->resinfo[resind].name);
177         }
178         else
179         {
180             rnm = *(atoms->resinfo[resind].rtp);
181         }
182
183         strcpy(atombuf, *(atoms->atomname[a]));
184         bReorderedNum = FALSE;
185         if (bReorderNum)
186         {
187             if (isdigit(atombuf[0]))
188             {
189                 c = atombuf[0];
190                 for (i = 0; ((size_t)i < strlen(atombuf)-1); i++)
191                 {
192                     atombuf[i] = atombuf[i+1];
193                 }
194                 atombuf[i]    = c;
195                 bReorderedNum = TRUE;
196             }
197         }
198         bRenamed = FALSE;
199         for (i = 0; (i < nxlate) && !bRenamed; i++)
200         {
201             /* Check if the base file name of the rtp and arn entry match */
202             if (restp == NULL ||
203                 gmx_strcasecmp(restp[resind].filebase, xlatom[i].filebase) == 0)
204             {
205                 /* Match the residue name */
206                 bMatch = (xlatom[i].res == NULL ||
207                           (gmx_strcasecmp("protein", xlatom[i].res) == 0 &&
208                            gmx_residuetype_is_protein(rt, rnm)) ||
209                           (gmx_strcasecmp("DNA", xlatom[i].res) == 0 &&
210                            gmx_residuetype_is_dna(rt, rnm)) ||
211                           (gmx_strcasecmp("RNA", xlatom[i].res) == 0 &&
212                            gmx_residuetype_is_rna(rt, rnm)));
213                 if (!bMatch)
214                 {
215                     ptr0 = rnm;
216                     ptr1 = xlatom[i].res;
217                     while (ptr0[0] != '\0' && ptr1[0] != '\0' &&
218                            (ptr0[0] == ptr1[0] || ptr1[0] == '?'))
219                     {
220                         ptr0++;
221                         ptr1++;
222                     }
223                     bMatch = (ptr0[0] == '\0' && ptr1[0] == '\0');
224                 }
225                 if (bMatch && strcmp(atombuf, xlatom[i].atom) == 0)
226                 {
227                     /* We have a match. */
228                     /* Don't free the old atomname,
229                      * since it might be in the symtab.
230                      */
231                     ptr0 = strdup(xlatom[i].replace);
232                     if (bVerbose)
233                     {
234                         printf("Renaming atom '%s' in residue %d %s to '%s'\n",
235                                *atoms->atomname[a],
236                                atoms->resinfo[resind].nr,
237                                *atoms->resinfo[resind].name,
238                                ptr0);
239                     }
240                     atoms->atomname[a] = put_symtab(symtab, ptr0);
241                     bRenamed           = TRUE;
242                 }
243             }
244         }
245         if (bReorderedNum && !bRenamed)
246         {
247             atoms->atomname[a] = put_symtab(symtab, atombuf);
248         }
249     }
250
251     done_xlatom(nxlate, xlatom);
252 }