2d8910747661e58305466d4cf1c47aa8ca862542
[alexxy/gromacs.git] / src / kernel / xlate.c
1 /*
2  * 
3  *                This source code is part of
4  * 
5  *                 G   R   O   M   A   C   S
6  * 
7  *          GROningen MAchine for Chemical Simulations
8  * 
9  *                        VERSION 3.2.0
10  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12  * Copyright (c) 2001-2004, The GROMACS development team,
13  * check out http://www.gromacs.org for more information.
14
15  * This program is free software; you can redistribute it and/or
16  * modify it under the terms of the GNU General Public License
17  * as published by the Free Software Foundation; either version 2
18  * of the License, or (at your option) any later version.
19  * 
20  * If you want to redistribute modifications, please consider that
21  * scientific software is very special. Version control is crucial -
22  * bugs must be traceable. We will be happy to consider code for
23  * inclusion in the official distribution, but derived work must not
24  * be called official GROMACS. Details are found in the README & COPYING
25  * files - if they are missing, get the official version at www.gromacs.org.
26  * 
27  * To help us fund GROMACS development, we humbly ask that you cite
28  * the papers on the package - you can find them in the top README file.
29  * 
30  * For more info, check our website at http://www.gromacs.org
31  * 
32  * And Hey:
33  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
34  */
35 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38
39 #include <ctype.h>
40 #include <string.h>
41 #include "typedefs.h"
42 #include "strdb.h"
43 #include "string2.h"
44 #include "smalloc.h"
45 #include "symtab.h"
46 #include "index.h"
47
48 typedef struct {
49   char *res;
50   char *atom;
51   char *replace;
52 } t_xlate_atom;
53
54 static t_xlate_atom *get_xlatoms(int *nxlatom)
55 {
56   const char  *xlfile="xlateat.dat";
57   
58   t_xlate_atom *xl=NULL;
59   char rbuf[32],abuf[32],repbuf[32];
60   char **lines,*_ptr;
61   int  nlines,i,n;
62   
63   nlines = get_lines(xlfile,&lines);
64   if (nlines > 0) 
65     snew(xl,nlines);
66     
67   n = 0;
68   for(i=0; (i<nlines); i++) {
69     if (sscanf(lines[i],"%s%s%s",rbuf,abuf,repbuf) != 3) 
70       fprintf(stderr,"Invalid line '%s' in %s\n",lines[i],xlfile);
71     else {
72       /* Use wildcards... */
73       if (strcmp(rbuf,"*") != 0)
74         xl[n].res = strdup(rbuf);
75       else
76         xl[n].res = NULL;
77       
78       /* Replace underscores in the string by spaces */
79       while ((_ptr = strchr(abuf,'_')) != 0)
80         *_ptr = ' ';
81       
82       xl[n].atom = strdup(abuf);
83       xl[n].replace = strdup(repbuf);
84       n++;
85     }
86     sfree(lines[i]);
87   }
88   if (nlines > 0)
89     sfree(lines);
90   fprintf(stderr,"%d out of %d lines of %s converted succesfully\n",
91           n,nlines,xlfile);
92   
93   *nxlatom = n;
94   
95   return xl;
96 }
97
98 static void done_xlatom(int nxlate,t_xlate_atom **xlatom)
99 {
100   int i;
101   
102   for(i=0; (i<nxlate); i++) {
103     if ((*xlatom)[i].res)
104       sfree((*xlatom)[i].res);
105     if ((*xlatom)[i].atom)
106       sfree((*xlatom)[i].atom);
107     if ((*xlatom)[i].replace)
108       sfree((*xlatom)[i].replace);
109   }
110   sfree(*xlatom);
111   *xlatom = NULL;
112 }
113
114 void rename_atoms(t_atoms *atoms,t_symtab *symtab,t_aa_names *aan)
115 {
116   int nxlate,a,i;
117   t_xlate_atom *xlatom;
118   char c,*res,atombuf[32];
119   bool bRenamed;
120
121   xlatom = get_xlatoms(&nxlate);
122
123   for(a=0; (a<atoms->nr); a++) {
124     res = *(atoms->resinfo[atoms->atom[a].resind].name);
125     strcpy(atombuf,*(atoms->atomname[a]));
126     if (isdigit(atombuf[0])) {
127       c = atombuf[0];
128       for (i=0; (i<strlen(atombuf)-1); i++)
129         atombuf[i]=atombuf[i+1];
130       atombuf[i]=c;
131     }
132     bRenamed=FALSE;
133     for(i=0; (i<nxlate) && !bRenamed; i++) {
134       if ((xlatom[i].res == NULL) || (strcasecmp(res,xlatom[i].res) == 0) ||
135           ((strcasecmp("protein",xlatom[i].res) == 0) && is_protein(aan,res)))
136         if (strcasecmp(atombuf,xlatom[i].atom) == 0) {
137           /* don't free the old atomname, since it might be in the symtab */
138           strcpy(atombuf,xlatom[i].replace);
139           bRenamed=TRUE;
140         }
141     }
142     if (strcmp(atombuf,*atoms->atomname[a]) != 0)
143       atoms->atomname[a] = put_symtab(symtab,atombuf);
144   }
145
146   done_xlatom(nxlate,&xlatom);
147 }
148