390f1dcef1bf517bae4b65b1e0d1986395625bc1
[alexxy/gromacs.git] / src / gromacs / topology / residuetypes.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2010,2013,2014,2015,2016 by the GROMACS development team.
5  * Copyright (c) 2017,2018,2019,2020,2021, 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.
9  *
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.
14  *
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.
19  *
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.
24  *
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.
32  *
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.
35  */
36 #include "gmxpre.h"
37
38 #include "residuetypes.h"
39
40 #include <cassert>
41 #include <cstdio>
42
43 #include <algorithm>
44 #include <optional>
45 #include <string>
46 #include <unordered_map>
47
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"
54
55 //! Definition for residue type that is not known.
56 const ResidueType c_undefinedResidueType = "Other";
57
58 //! Function object for comparisons used in std::unordered_map
59 class EqualCaseInsensitive
60 {
61 public:
62     bool operator()(const ResidueName& lhs, const ResidueName& rhs) const
63     {
64         return gmx::equalCaseInsensitive(lhs, rhs);
65     }
66 };
67
68 //! Implementation detail for ResidueTypeMap
69 class ResidueTypeMap::Impl
70 {
71 public:
72     //! Storage object for entries.
73     std::unordered_map<ResidueName, ResidueType, std::hash<ResidueName>, EqualCaseInsensitive> entries;
74 };
75
76 ResidueTypeMap::ResidueTypeMap() : impl_(new Impl)
77 {
78     char line[STRLEN];
79     char resname[STRLEN], restype[STRLEN], dum[STRLEN];
80
81     gmx::FilePtr db = gmx::openLibraryFile("residuetypes.dat");
82
83     while (get_a_line(db.get(), line, STRLEN))
84     {
85         strip_comment(line);
86         trim(line);
87         if (line[0] != '\0')
88         {
89             if (sscanf(line, "%1000s %1000s %1000s", resname, restype, dum) != 2)
90             {
91                 gmx_fatal(
92                         FARGS,
93                         "Incorrect number of columns (2 expected) for line in residuetypes.dat  ");
94             }
95             addResidue(resname, restype);
96         }
97     }
98 }
99
100 ResidueTypeMap::~ResidueTypeMap() = default;
101
102 bool ResidueTypeMap::nameIndexedInResidueTypeMap(const ResidueName& residueName)
103 {
104     return impl_->entries.find(residueName) != impl_->entries.end();
105 }
106
107 void ResidueTypeMap::addResidue(const ResidueName& residueName, const ResidueType& residueType)
108 {
109     if (auto [foundIt, insertionTookPlace] = impl_->entries.insert({ residueName, residueType });
110         !insertionTookPlace)
111     {
112         if (!gmx::equalCaseInsensitive(foundIt->second, residueType))
113         {
114             fprintf(stderr,
115                     "Warning: Residue '%s' already present with type '%s' in database, ignoring "
116                     "new type '%s'.\n",
117                     residueName.c_str(),
118                     foundIt->second.c_str(),
119                     residueType.c_str());
120         }
121     }
122 }
123
124 bool ResidueTypeMap::namedResidueHasType(const ResidueName& residueName, const ResidueType& residueType)
125 {
126     if (auto foundIt = impl_->entries.find(residueName); foundIt != impl_->entries.end())
127     {
128         return gmx::equalCaseInsensitive(residueType, foundIt->second);
129     }
130     return false;
131 }
132
133 ResidueType ResidueTypeMap::typeOfNamedDatabaseResidue(const ResidueName& residueName)
134 {
135     if (auto foundIt = impl_->entries.find(residueName); foundIt != impl_->entries.end())
136     {
137         return foundIt->second;
138     }
139     return c_undefinedResidueType;
140 }
141
142 std::optional<ResidueType> ResidueTypeMap::optionalTypeOfNamedDatabaseResidue(const ResidueName& residueName)
143 {
144     if (auto foundIt = impl_->entries.find(residueName); foundIt != impl_->entries.end())
145     {
146         return std::make_optional(foundIt->second);
147     }
148     return std::nullopt;
149 }