Code beautification with uncrustify
[alexxy/gromacs.git] / src / programs / pdb2gmx / xlate.c
1 /*  -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  *
4  *                This source code is part of
5  *
6  *                 G   R   O   M   A   C   S
7  *
8  *          GROningen MAchine for Chemical Simulations
9  *
10  *                        VERSION 3.2.0
11  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13  * Copyright (c) 2001-2004, The GROMACS development team,
14  * check out http://www.gromacs.org for more information.
15
16  * This program is free software; you can redistribute it and/or
17  * modify it under the terms of the GNU General Public License
18  * as published by the Free Software Foundation; either version 2
19  * of the License, or (at your option) any later version.
20  *
21  * If you want to redistribute modifications, please consider that
22  * scientific software is very special. Version control is crucial -
23  * bugs must be traceable. We will be happy to consider code for
24  * inclusion in the official distribution, but derived work must not
25  * be called official GROMACS. Details are found in the README & COPYING
26  * files - if they are missing, get the official version at www.gromacs.org.
27  *
28  * To help us fund GROMACS development, we humbly ask that you cite
29  * the papers on the package - you can find them in the top README file.
30  *
31  * For more info, check our website at http://www.gromacs.org
32  *
33  * And Hey:
34  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
35  */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include <ctype.h>
41 #include <string.h>
42 #include "typedefs.h"
43 #include "strdb.h"
44 #include "string2.h"
45 #include "smalloc.h"
46 #include "symtab.h"
47 #include "index.h"
48 #include "futil.h"
49 #include "fflibutil.h"
50 #include "hackblock.h"
51 #include "gmx_fatal.h"
52 #include "xlate.h"
53
54 typedef struct {
55     char *filebase;
56     char *res;
57     char *atom;
58     char *replace;
59 } t_xlate_atom;
60
61 static void get_xlatoms(const char *fn, FILE *fp,
62                         int *nptr, t_xlate_atom **xlptr)
63 {
64     char          filebase[STRLEN];
65     char          line[STRLEN];
66     char          abuf[1024], rbuf[1024], repbuf[1024], dumbuf[1024];
67     char         *_ptr;
68     int           n, na, idum;
69     t_xlate_atom *xl;
70
71     fflib_filename_base(fn, filebase, STRLEN);
72
73     n  = *nptr;
74     xl = *xlptr;
75
76     while (get_a_line(fp, line, STRLEN))
77     {
78         na = sscanf(line, "%s%s%s%s", rbuf, abuf, repbuf, dumbuf);
79         /* Check if we are reading an old format file with the number of items
80          * on the first line.
81          */
82         if (na == 1 && n == *nptr && sscanf(rbuf, "%d", &idum) == 1)
83         {
84             continue;
85         }
86         if (na != 3)
87         {
88             gmx_fatal(FARGS, "Expected a residue name and two atom names in file '%s', not '%s'", fn, line);
89         }
90
91         srenew(xl, n+1);
92         xl[n].filebase = strdup(filebase);
93
94         /* Use wildcards... */
95         if (strcmp(rbuf, "*") != 0)
96         {
97             xl[n].res = strdup(rbuf);
98         }
99         else
100         {
101             xl[n].res = NULL;
102         }
103
104         /* Replace underscores in the string by spaces */
105         while ((_ptr = strchr(abuf, '_')) != 0)
106         {
107             *_ptr = ' ';
108         }
109
110         xl[n].atom    = strdup(abuf);
111         xl[n].replace = strdup(repbuf);
112         n++;
113     }
114
115     *nptr  = n;
116     *xlptr = xl;
117 }
118
119 static void done_xlatom(int nxlate, t_xlate_atom *xlatom)
120 {
121     int i;
122
123     for (i = 0; (i < nxlate); i++)
124     {
125         sfree(xlatom[i].filebase);
126         if (xlatom[i].res != NULL)
127         {
128             sfree(xlatom[i].res);
129         }
130         sfree(xlatom[i].atom);
131         sfree(xlatom[i].replace);
132     }
133     sfree(xlatom);
134 }
135
136 void rename_atoms(const char *xlfile, const char *ffdir,
137                   t_atoms *atoms, t_symtab *symtab, const t_restp *restp,
138                   gmx_bool bResname, gmx_residuetype_t rt, gmx_bool bReorderNum,
139                   gmx_bool bVerbose)
140 {
141     FILE         *fp;
142     int           nxlate, a, i, resind;
143     t_xlate_atom *xlatom;
144     int           nf;
145     char        **f;
146     char          c, *rnm, atombuf[32], *ptr0, *ptr1;
147     gmx_bool      bReorderedNum, bRenamed, bMatch;
148
149     nxlate = 0;
150     xlatom = NULL;
151     if (xlfile != NULL)
152     {
153         fp = libopen(xlfile);
154         get_xlatoms(xlfile, fp, &nxlate, &xlatom);
155         fclose(fp);
156     }
157     else
158     {
159         nf = fflib_search_file_end(ffdir, ".arn", FALSE, &f);
160         for (i = 0; i < nf; i++)
161         {
162             fp = fflib_open(f[i]);
163             get_xlatoms(f[i], fp, &nxlate, &xlatom);
164             ffclose(fp);
165             sfree(f[i]);
166         }
167         sfree(f);
168     }
169
170     for (a = 0; (a < atoms->nr); a++)
171     {
172         resind = atoms->atom[a].resind;
173         if (bResname)
174         {
175             rnm = *(atoms->resinfo[resind].name);
176         }
177         else
178         {
179             rnm = *(atoms->resinfo[resind].rtp);
180         }
181
182         strcpy(atombuf, *(atoms->atomname[a]));
183         bReorderedNum = FALSE;
184         if (bReorderNum)
185         {
186             if (isdigit(atombuf[0]))
187             {
188                 c = atombuf[0];
189                 for (i = 0; ((size_t)i < strlen(atombuf)-1); i++)
190                 {
191                     atombuf[i] = atombuf[i+1];
192                 }
193                 atombuf[i]    = c;
194                 bReorderedNum = TRUE;
195             }
196         }
197         bRenamed = FALSE;
198         for (i = 0; (i < nxlate) && !bRenamed; i++)
199         {
200             /* Check if the base file name of the rtp and arn entry match */
201             if (restp == NULL ||
202                 gmx_strcasecmp(restp[resind].filebase, xlatom[i].filebase) == 0)
203             {
204                 /* Match the residue name */
205                 bMatch = (xlatom[i].res == NULL ||
206                           (gmx_strcasecmp("protein", xlatom[i].res) == 0 &&
207                            gmx_residuetype_is_protein(rt, rnm)) ||
208                           (gmx_strcasecmp("DNA", xlatom[i].res) == 0 &&
209                            gmx_residuetype_is_dna(rt, rnm)) ||
210                           (gmx_strcasecmp("RNA", xlatom[i].res) == 0 &&
211                            gmx_residuetype_is_rna(rt, rnm)));
212                 if (!bMatch)
213                 {
214                     ptr0 = rnm;
215                     ptr1 = xlatom[i].res;
216                     while (ptr0[0] != '\0' && ptr1[0] != '\0' &&
217                            (ptr0[0] == ptr1[0] || ptr1[0] == '?'))
218                     {
219                         ptr0++;
220                         ptr1++;
221                     }
222                     bMatch = (ptr0[0] == '\0' && ptr1[0] == '\0');
223                 }
224                 if (bMatch && strcmp(atombuf, xlatom[i].atom) == 0)
225                 {
226                     /* We have a match. */
227                     /* Don't free the old atomname,
228                      * since it might be in the symtab.
229                      */
230                     ptr0 = strdup(xlatom[i].replace);
231                     if (bVerbose)
232                     {
233                         printf("Renaming atom '%s' in residue %d %s to '%s'\n",
234                                *atoms->atomname[a],
235                                atoms->resinfo[resind].nr,
236                                *atoms->resinfo[resind].name,
237                                ptr0);
238                     }
239                     atoms->atomname[a] = put_symtab(symtab, ptr0);
240                     bRenamed           = TRUE;
241                 }
242             }
243         }
244         if (bReorderedNum && !bRenamed)
245         {
246             atoms->atomname[a] = put_symtab(symtab, atombuf);
247         }
248     }
249
250     done_xlatom(nxlate, xlatom);
251 }