Fixing copyright issues and code contributors
[alexxy/gromacs.git] / src / kernel / 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  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012,2013, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <ctype.h>
43 #include <string.h>
44 #include "typedefs.h"
45 #include "strdb.h"
46 #include "string2.h"
47 #include "smalloc.h"
48 #include "symtab.h"
49 #include "index.h"
50 #include "futil.h"
51 #include "fflibutil.h"
52 #include "hackblock.h"
53 #include "gmx_fatal.h"
54 #include "xlate.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             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             /* 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 }
253