Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / topology / residuetypes.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2010,2013,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 #include "gromacs/topology/residuetypes.h"
36
37 #include <cassert>
38 #include <cstdio>
39
40 #include "gromacs/fileio/strdb.h"
41 #include "gromacs/utility/cstringutil.h"
42 #include "gromacs/utility/fatalerror.h"
43 #include "gromacs/utility/futil.h"
44 #include "gromacs/utility/smalloc.h"
45
46 const char gmx_residuetype_undefined[] = "Other";
47
48 struct gmx_residuetype_t
49 {
50     int      n;
51     char **  resname;
52     char **  restype;
53 };
54
55 int
56 gmx_residuetype_init(gmx_residuetype_t **prt)
57 {
58     FILE                 *  db;
59     char                    line[STRLEN];
60     char                    resname[STRLEN], restype[STRLEN], dum[STRLEN];
61     gmx_residuetype_t      *rt;
62
63     snew(rt, 1);
64     *prt = rt;
65
66     rt->n        = 0;
67     rt->resname  = NULL;
68     rt->restype  = NULL;
69
70     db = libopen("residuetypes.dat");
71
72     while (get_a_line(db, line, STRLEN))
73     {
74         strip_comment(line);
75         trim(line);
76         if (line[0] != '\0')
77         {
78             if (sscanf(line, "%1000s %1000s %1000s", resname, restype, dum) != 2)
79             {
80                 gmx_fatal(FARGS, "Incorrect number of columns (2 expected) for line in residuetypes.dat");
81             }
82             gmx_residuetype_add(rt, resname, restype);
83         }
84     }
85
86     fclose(db);
87
88     return 0;
89 }
90
91 int
92 gmx_residuetype_destroy(gmx_residuetype_t *rt)
93 {
94     int i;
95
96     for (i = 0; i < rt->n; i++)
97     {
98         sfree(rt->resname[i]);
99         sfree(rt->restype[i]);
100     }
101     sfree(rt->resname);
102     sfree(rt->restype);
103     sfree(rt);
104
105     return 0;
106 }
107
108 /* Return 0 if the name was found, otherwise -1.
109  * p_restype is set to a pointer to the type name, or 'Other' if we did not find it.
110  */
111 int
112 gmx_residuetype_get_type(gmx_residuetype_t *rt, const char * resname, const char ** p_restype)
113 {
114     int    i, rc;
115
116     rc = -1;
117     for (i = 0; i < rt->n && rc; i++)
118     {
119         rc = gmx_strcasecmp(rt->resname[i], resname);
120     }
121
122     *p_restype = (rc == 0) ? rt->restype[i-1] : gmx_residuetype_undefined;
123
124     return rc;
125 }
126
127 int
128 gmx_residuetype_add(gmx_residuetype_t *rt, const char *newresname, const char *newrestype)
129 {
130     int           found;
131     const char *  p_oldtype;
132
133     found = !gmx_residuetype_get_type(rt, newresname, &p_oldtype);
134
135     if (found && gmx_strcasecmp(p_oldtype, newrestype))
136     {
137         fprintf(stderr, "Warning: Residue '%s' already present with type '%s' in database, ignoring new type '%s'.",
138                 newresname, p_oldtype, newrestype);
139     }
140
141     if (found == 0)
142     {
143         srenew(rt->resname, rt->n+1);
144         srenew(rt->restype, rt->n+1);
145         rt->resname[rt->n] = gmx_strdup(newresname);
146         rt->restype[rt->n] = gmx_strdup(newrestype);
147         rt->n++;
148     }
149
150     return 0;
151 }
152
153 int
154 gmx_residuetype_get_alltypes(gmx_residuetype_t   *rt,
155                              const char ***       p_typenames,
156                              int *                ntypes)
157 {
158     int            i, n;
159     const char **  my_typename;
160
161     n           = 0;
162     my_typename = NULL;
163     for (i = 0; i < rt->n; i++)
164     {
165         const char *const p      = rt->restype[i];
166         bool              bFound = false;
167         for (int j = 0; j < n && !bFound; j++)
168         {
169             assert(my_typename != NULL);
170             bFound = !gmx_strcasecmp(p, my_typename[j]);
171         }
172         if (!bFound)
173         {
174             srenew(my_typename, n+1);
175             my_typename[n] = p;
176             n++;
177         }
178     }
179     *ntypes      = n;
180     *p_typenames = my_typename;
181
182     return 0;
183 }
184
185 gmx_bool
186 gmx_residuetype_is_protein(gmx_residuetype_t *rt, const char *resnm)
187 {
188     gmx_bool    rc;
189     const char *p_type;
190
191     if (gmx_residuetype_get_type(rt, resnm, &p_type) == 0 &&
192         gmx_strcasecmp(p_type, "Protein") == 0)
193     {
194         rc = TRUE;
195     }
196     else
197     {
198         rc = FALSE;
199     }
200     return rc;
201 }
202
203 gmx_bool
204 gmx_residuetype_is_dna(gmx_residuetype_t *rt, const char *resnm)
205 {
206     gmx_bool    rc;
207     const char *p_type;
208
209     if (gmx_residuetype_get_type(rt, resnm, &p_type) == 0 &&
210         gmx_strcasecmp(p_type, "DNA") == 0)
211     {
212         rc = TRUE;
213     }
214     else
215     {
216         rc = FALSE;
217     }
218     return rc;
219 }
220
221 gmx_bool
222 gmx_residuetype_is_rna(gmx_residuetype_t *rt, const char *resnm)
223 {
224     gmx_bool    rc;
225     const char *p_type;
226
227     if (gmx_residuetype_get_type(rt, resnm, &p_type) == 0 &&
228         gmx_strcasecmp(p_type, "RNA") == 0)
229     {
230         rc = TRUE;
231     }
232     else
233     {
234         rc = FALSE;
235     }
236     return rc;
237 }
238
239 /* Return the size of the arrays */
240 int
241 gmx_residuetype_get_size(gmx_residuetype_t *rt)
242 {
243     return rt->n;
244 }
245
246 /* Search for a residuetype with name resnm within the
247  * gmx_residuetype database. Return the index if found,
248  * otherwise -1.
249  */
250 int
251 gmx_residuetype_get_index(gmx_residuetype_t *rt, const char *resnm)
252 {
253     int i, rc;
254
255     rc = -1;
256     for (i = 0; i < rt->n && rc; i++)
257     {
258         rc = gmx_strcasecmp(rt->resname[i], resnm);
259     }
260
261     return (0 == rc) ? i-1 : -1;
262 }
263
264 /* Return the name of the residuetype with the given index, or
265  * NULL if not found. */
266 const char *
267 gmx_residuetype_get_name(gmx_residuetype_t *rt, int index)
268 {
269     if (index >= 0 && index < rt->n)
270     {
271         return rt->resname[index];
272     }
273     else
274     {
275         return NULL;
276     }
277 }