Update instructions in containers.rst
[gromacs.git] / src / gromacs / topology / residuetypes.cpp
blob73aad6d57f7f7dd9bf4b0b4f27c2620eef0af61c
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.
5 * Copyright (c) 2017,2018,2019,2020, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
36 #include "gmxpre.h"
38 #include "residuetypes.h"
40 #include <cassert>
41 #include <cstdio>
43 #include <algorithm>
44 #include <iterator>
45 #include <optional>
46 #include <string>
48 #include "gromacs/utility/arrayref.h"
49 #include "gromacs/utility/cstringutil.h"
50 #include "gromacs/utility/fatalerror.h"
51 #include "gromacs/utility/futil.h"
52 #include "gromacs/utility/smalloc.h"
53 #include "gromacs/utility/strdb.h"
55 //! Definition for residue type that is not known.
56 const std::string c_undefinedResidueType = "Other";
58 //! Single ResidueType entry object
59 struct ResidueTypeEntry
61 //! Default constructor creates complete object.
62 ResidueTypeEntry(const std::string& rName, const std::string& rType) :
63 residueName(rName),
64 residueType(rType)
67 //! Name of the residue in the entry.
68 std::string residueName;
69 //! Type of the residue in the entry.
70 std::string residueType;
73 //! Implementation detail for ResidueTypes
74 class ResidueType::Impl
76 public:
77 //! Storage object for entries.
78 std::vector<ResidueTypeEntry> entry;
81 ResidueType::ResidueType() : impl_(new Impl)
83 char line[STRLEN];
84 char resname[STRLEN], restype[STRLEN], dum[STRLEN];
86 gmx::FilePtr db = gmx::openLibraryFile("residuetypes.dat");
88 while (get_a_line(db.get(), line, STRLEN))
90 strip_comment(line);
91 trim(line);
92 if (line[0] != '\0')
94 if (sscanf(line, "%1000s %1000s %1000s", resname, restype, dum) != 2)
96 gmx_fatal(
97 FARGS,
98 "Incorrect number of columns (2 expected) for line in residuetypes.dat ");
100 addResidue(resname, restype);
105 ResidueType::~ResidueType() {}
107 /*! \brief
108 * Return an optional const iterator to a residue entry that matches the given name.
110 * \param[in] entries Currently registered residue entries in the database.
111 * \param[in] residueName Name of a residue to compare to database.
112 * \returns An optional iterator to the residue entry that was found.
114 static std::optional<gmx::ArrayRef<const ResidueTypeEntry>::const_iterator>
115 findResidueEntryWithName(gmx::ArrayRef<const ResidueTypeEntry> entries, const std::string& residueName)
117 auto foundIt =
118 std::find_if(entries.begin(), entries.end(), [&residueName](const ResidueTypeEntry& old) {
119 return gmx::equalCaseInsensitive(residueName, old.residueName);
121 return (foundIt != entries.end()) ? std::make_optional(foundIt) : std::nullopt;
124 bool ResidueType::nameIndexedInResidueTypes(const std::string& residueName)
126 return findResidueEntryWithName(impl_->entry, residueName).has_value();
129 void ResidueType::addResidue(const std::string& residueName, const std::string& residueType)
131 if (auto foundIt = findResidueEntryWithName(impl_->entry, residueName))
133 if (!gmx::equalCaseInsensitive((*foundIt)->residueType, residueType))
135 fprintf(stderr,
136 "Warning: Residue '%s' already present with type '%s' in database, ignoring "
137 "new type '%s'.\n",
138 residueName.c_str(), (*foundIt)->residueType.c_str(), residueType.c_str());
141 else
143 impl_->entry.emplace_back(residueName, residueType);
147 bool ResidueType::namedResidueHasType(const std::string& residueName, const std::string& residueType)
149 auto foundIt = findResidueEntryWithName(impl_->entry, residueName);
150 return foundIt ? gmx::equalCaseInsensitive(residueType, (*foundIt)->residueType) : false;
153 int ResidueType::numberOfEntries() const
155 return impl_->entry.size();
158 int ResidueType::indexFromResidueName(const std::string& residueName) const
160 gmx::ArrayRef<const ResidueTypeEntry> temp(impl_->entry);
161 auto foundIt = findResidueEntryWithName(temp, residueName);
162 return foundIt ? std::distance(temp.begin(), *foundIt) : -1;
165 std::string ResidueType::nameFromResidueIndex(int index) const
167 if (index >= 0 && index < gmx::ssize(impl_->entry))
169 return impl_->entry[index].residueName;
171 else
173 return "";
177 std::string ResidueType::typeOfNamedDatabaseResidue(const std::string& residueName)
179 auto foundIt = findResidueEntryWithName(impl_->entry, residueName);
180 return foundIt ? (*foundIt)->residueType : c_undefinedResidueType;
183 std::optional<std::string> ResidueType::optionalTypeOfNamedDatabaseResidue(const std::string& residueName)
185 auto foundIt = findResidueEntryWithName(impl_->entry, residueName);
186 return foundIt ? std::make_optional((*foundIt)->residueType) : std::nullopt;