Added pull coordinate geometry angle-axis
[gromacs.git] / src / gromacs / topology / residuetypes.cpp
blob50485b318f8c0c1365945114c4019d39802fbdb9
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2010,2013,2014,2015,2016, 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.
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.
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.
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.
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.
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.
35 #include "gmxpre.h"
37 #include "residuetypes.h"
39 #include <cassert>
40 #include <cstdio>
42 #include "gromacs/utility/cstringutil.h"
43 #include "gromacs/utility/fatalerror.h"
44 #include "gromacs/utility/futil.h"
45 #include "gromacs/utility/smalloc.h"
46 #include "gromacs/utility/strdb.h"
48 const char gmx_residuetype_undefined[] = "Other";
50 struct gmx_residuetype_t
52 int n;
53 char ** resname;
54 char ** restype;
57 int
58 gmx_residuetype_init(gmx_residuetype_t **prt)
60 FILE * db;
61 char line[STRLEN];
62 char resname[STRLEN], restype[STRLEN], dum[STRLEN];
63 gmx_residuetype_t *rt;
65 snew(rt, 1);
66 *prt = rt;
68 rt->n = 0;
69 rt->resname = NULL;
70 rt->restype = NULL;
72 db = libopen("residuetypes.dat");
74 while (get_a_line(db, line, STRLEN))
76 strip_comment(line);
77 trim(line);
78 if (line[0] != '\0')
80 if (sscanf(line, "%1000s %1000s %1000s", resname, restype, dum) != 2)
82 gmx_fatal(FARGS, "Incorrect number of columns (2 expected) for line in residuetypes.dat");
84 gmx_residuetype_add(rt, resname, restype);
88 fclose(db);
90 return 0;
93 int
94 gmx_residuetype_destroy(gmx_residuetype_t *rt)
96 int i;
98 for (i = 0; i < rt->n; i++)
100 sfree(rt->resname[i]);
101 sfree(rt->restype[i]);
103 sfree(rt->resname);
104 sfree(rt->restype);
105 sfree(rt);
107 return 0;
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.
114 gmx_residuetype_get_type(gmx_residuetype_t *rt, const char * resname, const char ** p_restype)
116 int i, rc;
118 rc = -1;
119 for (i = 0; i < rt->n && rc; i++)
121 rc = gmx_strcasecmp(rt->resname[i], resname);
124 *p_restype = (rc == 0) ? rt->restype[i-1] : gmx_residuetype_undefined;
126 return rc;
130 gmx_residuetype_add(gmx_residuetype_t *rt, const char *newresname, const char *newrestype)
132 int found;
133 const char * p_oldtype;
135 found = !gmx_residuetype_get_type(rt, newresname, &p_oldtype);
137 if (found && gmx_strcasecmp(p_oldtype, newrestype))
139 fprintf(stderr, "Warning: Residue '%s' already present with type '%s' in database, ignoring new type '%s'.",
140 newresname, p_oldtype, newrestype);
143 if (found == 0)
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++;
152 return 0;
156 gmx_residuetype_get_alltypes(gmx_residuetype_t *rt,
157 const char *** p_typenames,
158 int * ntypes)
160 int n = 0;
161 const char **my_typename = NULL;
163 if (rt->n > 0)
165 int i = 0;
166 const char *p = rt->restype[i];
167 snew(my_typename, n+1);
168 my_typename[n] = p;
169 n = 1;
171 for (i = 1; i < rt->n; i++)
173 p = rt->restype[i];
174 bool bFound = false;
175 for (int j = 0; j < n && !bFound; j++)
177 bFound = !gmx_strcasecmp(p, my_typename[j]);
179 if (!bFound)
181 srenew(my_typename, n+1);
182 my_typename[n] = p;
183 n++;
187 *ntypes = n;
188 *p_typenames = my_typename;
190 return 0;
193 gmx_bool
194 gmx_residuetype_is_protein(gmx_residuetype_t *rt, const char *resnm)
196 gmx_bool rc;
197 const char *p_type;
199 if (gmx_residuetype_get_type(rt, resnm, &p_type) == 0 &&
200 gmx_strcasecmp(p_type, "Protein") == 0)
202 rc = TRUE;
204 else
206 rc = FALSE;
208 return rc;
211 gmx_bool
212 gmx_residuetype_is_dna(gmx_residuetype_t *rt, const char *resnm)
214 gmx_bool rc;
215 const char *p_type;
217 if (gmx_residuetype_get_type(rt, resnm, &p_type) == 0 &&
218 gmx_strcasecmp(p_type, "DNA") == 0)
220 rc = TRUE;
222 else
224 rc = FALSE;
226 return rc;
229 gmx_bool
230 gmx_residuetype_is_rna(gmx_residuetype_t *rt, const char *resnm)
232 gmx_bool rc;
233 const char *p_type;
235 if (gmx_residuetype_get_type(rt, resnm, &p_type) == 0 &&
236 gmx_strcasecmp(p_type, "RNA") == 0)
238 rc = TRUE;
240 else
242 rc = FALSE;
244 return rc;
247 /* Return the size of the arrays */
249 gmx_residuetype_get_size(gmx_residuetype_t *rt)
251 return rt->n;
254 /* Search for a residuetype with name resnm within the
255 * gmx_residuetype database. Return the index if found,
256 * otherwise -1.
259 gmx_residuetype_get_index(gmx_residuetype_t *rt, const char *resnm)
261 int i, rc;
263 rc = -1;
264 for (i = 0; i < rt->n && rc; i++)
266 rc = gmx_strcasecmp(rt->resname[i], resnm);
269 return (0 == rc) ? i-1 : -1;
272 /* Return the name of the residuetype with the given index, or
273 * NULL if not found. */
274 const char *
275 gmx_residuetype_get_name(gmx_residuetype_t *rt, int index)
277 if (index >= 0 && index < rt->n)
279 return rt->resname[index];
281 else
283 return NULL;