Sort all includes in src/gromacs
[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 "gmxpre.h"
36
37 #include "residuetypes.h"
38
39 #include <cassert>
40 #include <cstdio>
41
42 #include "gromacs/fileio/strdb.h"
43 #include "gromacs/utility/cstringutil.h"
44 #include "gromacs/utility/fatalerror.h"
45 #include "gromacs/utility/futil.h"
46 #include "gromacs/utility/smalloc.h"
47
48 const char gmx_residuetype_undefined[] = "Other";
49
50 struct gmx_residuetype_t
51 {
52     int      n;
53     char **  resname;
54     char **  restype;
55 };
56
57 int
58 gmx_residuetype_init(gmx_residuetype_t **prt)
59 {
60     FILE                 *  db;
61     char                    line[STRLEN];
62     char                    resname[STRLEN], restype[STRLEN], dum[STRLEN];
63     gmx_residuetype_t      *rt;
64
65     snew(rt, 1);
66     *prt = rt;
67
68     rt->n        = 0;
69     rt->resname  = NULL;
70     rt->restype  = NULL;
71
72     db = libopen("residuetypes.dat");
73
74     while (get_a_line(db, line, STRLEN))
75     {
76         strip_comment(line);
77         trim(line);
78         if (line[0] != '\0')
79         {
80             if (sscanf(line, "%1000s %1000s %1000s", resname, restype, dum) != 2)
81             {
82                 gmx_fatal(FARGS, "Incorrect number of columns (2 expected) for line in residuetypes.dat");
83             }
84             gmx_residuetype_add(rt, resname, restype);
85         }
86     }
87
88     fclose(db);
89
90     return 0;
91 }
92
93 int
94 gmx_residuetype_destroy(gmx_residuetype_t *rt)
95 {
96     int i;
97
98     for (i = 0; i < rt->n; i++)
99     {
100         sfree(rt->resname[i]);
101         sfree(rt->restype[i]);
102     }
103     sfree(rt->resname);
104     sfree(rt->restype);
105     sfree(rt);
106
107     return 0;
108 }
109
110 /* Return 0 if the name was found, otherwise -1.
111  * p_restype is set to a pointer to the type name, or 'Other' if we did not find it.
112  */
113 int
114 gmx_residuetype_get_type(gmx_residuetype_t *rt, const char * resname, const char ** p_restype)
115 {
116     int    i, rc;
117
118     rc = -1;
119     for (i = 0; i < rt->n && rc; i++)
120     {
121         rc = gmx_strcasecmp(rt->resname[i], resname);
122     }
123
124     *p_restype = (rc == 0) ? rt->restype[i-1] : gmx_residuetype_undefined;
125
126     return rc;
127 }
128
129 int
130 gmx_residuetype_add(gmx_residuetype_t *rt, const char *newresname, const char *newrestype)
131 {
132     int           found;
133     const char *  p_oldtype;
134
135     found = !gmx_residuetype_get_type(rt, newresname, &p_oldtype);
136
137     if (found && gmx_strcasecmp(p_oldtype, newrestype))
138     {
139         fprintf(stderr, "Warning: Residue '%s' already present with type '%s' in database, ignoring new type '%s'.",
140                 newresname, p_oldtype, newrestype);
141     }
142
143     if (found == 0)
144     {
145         srenew(rt->resname, rt->n+1);
146         srenew(rt->restype, rt->n+1);
147         rt->resname[rt->n] = gmx_strdup(newresname);
148         rt->restype[rt->n] = gmx_strdup(newrestype);
149         rt->n++;
150     }
151
152     return 0;
153 }
154
155 int
156 gmx_residuetype_get_alltypes(gmx_residuetype_t   *rt,
157                              const char ***       p_typenames,
158                              int *                ntypes)
159 {
160     int            i, n;
161     const char **  my_typename;
162
163     n           = 0;
164     my_typename = NULL;
165     for (i = 0; i < rt->n; i++)
166     {
167         const char *const p      = rt->restype[i];
168         bool              bFound = false;
169         for (int j = 0; j < n && !bFound; j++)
170         {
171             assert(my_typename != NULL);
172             bFound = !gmx_strcasecmp(p, my_typename[j]);
173         }
174         if (!bFound)
175         {
176             srenew(my_typename, n+1);
177             my_typename[n] = p;
178             n++;
179         }
180     }
181     *ntypes      = n;
182     *p_typenames = my_typename;
183
184     return 0;
185 }
186
187 gmx_bool
188 gmx_residuetype_is_protein(gmx_residuetype_t *rt, const char *resnm)
189 {
190     gmx_bool    rc;
191     const char *p_type;
192
193     if (gmx_residuetype_get_type(rt, resnm, &p_type) == 0 &&
194         gmx_strcasecmp(p_type, "Protein") == 0)
195     {
196         rc = TRUE;
197     }
198     else
199     {
200         rc = FALSE;
201     }
202     return rc;
203 }
204
205 gmx_bool
206 gmx_residuetype_is_dna(gmx_residuetype_t *rt, const char *resnm)
207 {
208     gmx_bool    rc;
209     const char *p_type;
210
211     if (gmx_residuetype_get_type(rt, resnm, &p_type) == 0 &&
212         gmx_strcasecmp(p_type, "DNA") == 0)
213     {
214         rc = TRUE;
215     }
216     else
217     {
218         rc = FALSE;
219     }
220     return rc;
221 }
222
223 gmx_bool
224 gmx_residuetype_is_rna(gmx_residuetype_t *rt, const char *resnm)
225 {
226     gmx_bool    rc;
227     const char *p_type;
228
229     if (gmx_residuetype_get_type(rt, resnm, &p_type) == 0 &&
230         gmx_strcasecmp(p_type, "RNA") == 0)
231     {
232         rc = TRUE;
233     }
234     else
235     {
236         rc = FALSE;
237     }
238     return rc;
239 }
240
241 /* Return the size of the arrays */
242 int
243 gmx_residuetype_get_size(gmx_residuetype_t *rt)
244 {
245     return rt->n;
246 }
247
248 /* Search for a residuetype with name resnm within the
249  * gmx_residuetype database. Return the index if found,
250  * otherwise -1.
251  */
252 int
253 gmx_residuetype_get_index(gmx_residuetype_t *rt, const char *resnm)
254 {
255     int i, rc;
256
257     rc = -1;
258     for (i = 0; i < rt->n && rc; i++)
259     {
260         rc = gmx_strcasecmp(rt->resname[i], resnm);
261     }
262
263     return (0 == rc) ? i-1 : -1;
264 }
265
266 /* Return the name of the residuetype with the given index, or
267  * NULL if not found. */
268 const char *
269 gmx_residuetype_get_name(gmx_residuetype_t *rt, int index)
270 {
271     if (index >= 0 && index < rt->n)
272     {
273         return rt->resname[index];
274     }
275     else
276     {
277         return NULL;
278     }
279 }