Add gmx convert-trj
[gromacs.git] / src / gromacs / topology / residuetypes.cpp
blob8387f8fe33ffae74e24adb3059fcc6761b51e589
1 /*
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.
35 #include "gmxpre.h"
37 #include "residuetypes.h"
39 #include <cassert>
40 #include <cstdio>
42 #include <algorithm>
43 #include <iterator>
44 #include <string>
46 #include "gromacs/compat/optional.h"
47 #include "gromacs/utility/arrayref.h"
48 #include "gromacs/utility/cstringutil.h"
49 #include "gromacs/utility/fatalerror.h"
50 #include "gromacs/utility/futil.h"
51 #include "gromacs/utility/smalloc.h"
52 #include "gromacs/utility/strdb.h"
54 //! Definition for residue type that is not known.
55 const std::string c_undefinedResidueType = "Other";
57 //! Single ResidueType entry object
58 struct ResidueTypeEntry
60 //! Default constructor creates complete object.
61 ResidueTypeEntry(const std::string &rName, const std::string &rType)
62 : residueName(rName), residueType(rType)
64 //! Name of the residue in the entry.
65 std::string residueName;
66 //! Type of the residue in the entry.
67 std::string residueType;
70 //! Implementation detail for ResidueTypes
71 class ResidueType::Impl
73 public:
74 //! Storage object for entries.
75 std::vector<ResidueTypeEntry> entry;
78 ResidueType::ResidueType()
79 : impl_(new Impl)
81 char line[STRLEN];
82 char resname[STRLEN], restype[STRLEN], dum[STRLEN];
84 gmx::FilePtr db = gmx::openLibraryFile("residuetypes.dat");
86 while (get_a_line(db.get(), line, STRLEN))
88 strip_comment(line);
89 trim(line);
90 if (line[0] != '\0')
92 if (sscanf(line, "%1000s %1000s %1000s", resname, restype, dum) != 2)
94 gmx_fatal(FARGS, "Incorrect number of columns (2 expected) for line in residuetypes.dat ");
96 addResidue(resname, restype);
101 ResidueType::~ResidueType()
105 /*! \brief
106 * Return an optional const iterator to a residue entry that matches the given name.
108 * \param[in] entries Currently registered residue entries in the database.
109 * \param[in] residueName Name of a residue to compare to database.
110 * \returns An optional iterator to the residue entry that was found.
112 static gmx::compat::optional<gmx::ArrayRef<const ResidueTypeEntry>::const_iterator>
113 findResidueEntryWithName(gmx::ArrayRef<const ResidueTypeEntry> entries, const std::string &residueName)
115 auto foundIt = std::find_if(entries.begin(), entries.end(),
116 [&residueName](const ResidueTypeEntry &old)
117 { return gmx::equalCaseInsensitive(residueName, old.residueName); });
118 return (foundIt != entries.end()) ? gmx::compat::make_optional(foundIt) : gmx::compat::nullopt;
121 bool ResidueType::nameIndexedInResidueTypes(const std::string &residueName)
123 return findResidueEntryWithName(impl_->entry, residueName).has_value();
126 void ResidueType::addResidue(const std::string &residueName, const std::string &residueType)
128 if (auto foundIt = findResidueEntryWithName(impl_->entry, residueName))
130 if (!gmx::equalCaseInsensitive((*foundIt)->residueType, residueType))
132 fprintf(stderr, "Warning: Residue '%s' already present with type '%s' in database, ignoring new type '%s'.\n",
133 residueName.c_str(), (*foundIt)->residueType.c_str(), residueType.c_str());
136 else
138 impl_->entry.emplace_back(residueName, residueType);
142 bool ResidueType::namedResidueHasType(const std::string &residueName, const std::string &residueType)
144 auto foundIt = findResidueEntryWithName(impl_->entry, residueName);
145 return foundIt ? gmx::equalCaseInsensitive(residueType, (*foundIt)->residueType) : false;
148 int ResidueType::numberOfEntries() const
150 return impl_->entry.size();
153 int ResidueType::indexFromResidueName(const std::string &residueName) const
155 gmx::ArrayRef<const ResidueTypeEntry> temp(impl_->entry);
156 auto foundIt = findResidueEntryWithName(temp, residueName);
157 return foundIt ? std::distance(temp.begin(), *foundIt) : -1;
160 std::string ResidueType::nameFromResidueIndex(int index) const
162 if (index >= 0 && index < gmx::ssize(impl_->entry))
164 return impl_->entry[index].residueName;
166 else
168 return "";
172 std::string
173 ResidueType::typeOfNamedDatabaseResidue(const std::string &residueName)
175 auto foundIt = findResidueEntryWithName(impl_->entry, residueName);
176 return foundIt ? (*foundIt)->residueType : c_undefinedResidueType;
179 gmx::compat::optional<std::string>
180 ResidueType::optionalTypeOfNamedDatabaseResidue(const std::string &residueName)
182 auto foundIt = findResidueEntryWithName(impl_->entry, residueName);
183 return foundIt ? gmx::compat::make_optional((*foundIt)->residueType) : gmx::compat::nullopt;