2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2010,2013,2014,2015,2016,2017, 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.
37 #include "residuetypes.h"
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
58 gmx_residuetype_init(gmx_residuetype_t
**prt
)
62 char resname
[STRLEN
], restype
[STRLEN
], dum
[STRLEN
];
63 gmx_residuetype_t
*rt
;
69 rt
->resname
= nullptr;
70 rt
->restype
= nullptr;
72 db
= libopen("residuetypes.dat");
74 while (get_a_line(db
, line
, STRLEN
))
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
);
94 gmx_residuetype_destroy(gmx_residuetype_t
*rt
)
98 for (i
= 0; i
< rt
->n
; i
++)
100 sfree(rt
->resname
[i
]);
101 sfree(rt
->restype
[i
]);
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
)
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
;
130 gmx_residuetype_add(gmx_residuetype_t
*rt
, const char *newresname
, const char *newrestype
)
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
);
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
);
156 gmx_residuetype_get_alltypes(gmx_residuetype_t
*rt
,
157 const char *** p_typenames
,
161 const char **my_typename
= nullptr;
166 const char *p
= rt
->restype
[i
];
167 snew(my_typename
, n
+1);
171 for (i
= 1; i
< rt
->n
; i
++)
175 for (int j
= 0; j
< n
&& !bFound
; j
++)
177 bFound
= !gmx_strcasecmp(p
, my_typename
[j
]);
181 srenew(my_typename
, n
+1);
188 *p_typenames
= my_typename
;
194 gmx_residuetype_is_protein(gmx_residuetype_t
*rt
, const char *resnm
)
199 if (gmx_residuetype_get_type(rt
, resnm
, &p_type
) == 0 &&
200 gmx_strcasecmp(p_type
, "Protein") == 0)
212 gmx_residuetype_is_dna(gmx_residuetype_t
*rt
, const char *resnm
)
217 if (gmx_residuetype_get_type(rt
, resnm
, &p_type
) == 0 &&
218 gmx_strcasecmp(p_type
, "DNA") == 0)
230 gmx_residuetype_is_rna(gmx_residuetype_t
*rt
, const char *resnm
)
235 if (gmx_residuetype_get_type(rt
, resnm
, &p_type
) == 0 &&
236 gmx_strcasecmp(p_type
, "RNA") == 0)
247 /* Return the size of the arrays */
249 gmx_residuetype_get_size(gmx_residuetype_t
*rt
)
254 /* Search for a residuetype with name resnm within the
255 * gmx_residuetype database. Return the index if found,
259 gmx_residuetype_get_index(gmx_residuetype_t
*rt
, const char *resnm
)
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. */
275 gmx_residuetype_get_name(gmx_residuetype_t
*rt
, int index
)
277 if (index
>= 0 && index
< rt
->n
)
279 return rt
->resname
[index
];