2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2010,2013,2014,2015,2016,2017,2018,2019, 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"
46 #include "gromacs/utility/arrayref.h"
47 #include "gromacs/utility/cstringutil.h"
48 #include "gromacs/utility/fatalerror.h"
49 #include "gromacs/utility/futil.h"
50 #include "gromacs/utility/smalloc.h"
51 #include "gromacs/utility/strdb.h"
53 //! Definition for residue type that is not known.
54 const std::string c_undefinedResidueType = "Other";
56 //! Single ResidueType entry object
57 struct ResidueTypeEntry
59 //! Default constructor creates complete object.
60 ResidueTypeEntry(const std::string &rName, const std::string &rType)
61 : residueName(rName), residueType(rType)
63 //! Name of the residue in the entry.
64 std::string residueName;
65 //! Type of the residue in the entry.
66 std::string residueType;
69 //! Implementation detail for ResidueTypes
70 class ResidueType::Impl
73 //! Storage object for entries.
74 std::vector<ResidueTypeEntry> entry;
77 ResidueType::ResidueType()
81 char resname[STRLEN], restype[STRLEN], dum[STRLEN];
83 gmx::FilePtr db = gmx::openLibraryFile("residuetypes.dat");
85 while (get_a_line(db.get(), line, STRLEN))
91 if (sscanf(line, "%1000s %1000s %1000s", resname, restype, dum) != 2)
93 gmx_fatal(FARGS, "Incorrect number of columns (2 expected) for line in residuetypes.dat ");
95 addResidue(resname, restype);
100 ResidueType::~ResidueType()
105 * Find a residue entry by the residue name.
107 * \param[in] entries Currently registered residue entries in the database.
108 * \param[in] residueName Name of a residue to compare to database.
109 * \returns A pointer to the entry that was found, or nullptr.
111 static gmx::ArrayRef<const ResidueTypeEntry>::iterator
112 residueEntryByResidueName(gmx::ArrayRef<const ResidueTypeEntry> entries, const std::string &residueName)
114 return std::find_if(entries.begin(), entries.end(),
115 [&residueName](const ResidueTypeEntry &old)
116 { return gmx::equalCaseInsensitive(residueName, old.residueName); });
119 bool ResidueType::nameIndexedInResidueTypes(const std::string &residueName)
121 return residueEntryByResidueName(impl_->entry, residueName) != nullptr;
124 void ResidueType::addResidue(const std::string &residueName, const std::string &residueType)
126 gmx::ArrayRef<const ResidueTypeEntry> temp(impl_->entry);
127 auto found = residueEntryByResidueName(temp, residueName);
129 if (found != temp.end())
131 if (!gmx::equalCaseInsensitive(found->residueType, residueType))
133 fprintf(stderr, "Warning: Residue '%s' already present with type '%s' in database, ignoring new type '%s'.\n",
134 residueName.c_str(), found->residueType.c_str(), residueType.c_str());
139 impl_->entry.emplace_back(residueName, residueType);
143 bool ResidueType::namedResidueHasType(const std::string &residueName, const std::string &residueType)
145 gmx::ArrayRef<const ResidueTypeEntry> temp(impl_->entry);
146 auto found = residueEntryByResidueName(temp, residueName);
147 return ((found != temp.end()) &&
148 gmx::equalCaseInsensitive(residueType, found->residueType));
151 int ResidueType::numberOfEntries() const
153 return impl_->entry.size();
156 int ResidueType::indexFromResidueName(const std::string &residueName) const
158 gmx::ArrayRef<const ResidueTypeEntry> temp(impl_->entry);
159 auto found = residueEntryByResidueName(temp, residueName);
160 return (found != temp.end()) ? std::distance(temp.begin(), found) : -1;
163 const std::string ResidueType::nameFromResidueIndex(int index) const
165 if (index >= 0 && index < gmx::ssize(impl_->entry))
167 return impl_->entry[index].residueName;
175 const std::string ResidueType::typeNameForIndexedResidue(const std::string &residueName)
177 gmx::ArrayRef<const ResidueTypeEntry> temp(impl_->entry);
178 auto found = residueEntryByResidueName(temp, residueName);
179 if (found != temp.end())
181 return found->residueType;
185 return c_undefinedResidueType;