Add function to get residue start and end
[gromacs.git] / src / gromacs / selection / symrec.cpp
blobb435dda6b7ea1178c125c2dfeaaaf44b244ca8af
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009-2017, The GROMACS development team.
5 * Copyright (c) 2019, 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 /*! \internal \file
37 * \brief
38 * Implements classes in symrec.h.
40 * \author Teemu Murtola <teemu.murtola@gmail.com>
41 * \ingroup module_selection
43 #include "gmxpre.h"
45 #include "symrec.h"
47 #include <map>
48 #include <memory>
49 #include <string>
50 #include <utility>
52 #include "gromacs/utility/arraysize.h"
53 #include "gromacs/utility/exceptions.h"
54 #include "gromacs/utility/gmxassert.h"
55 #include "gromacs/utility/stringutil.h"
57 #include "poscalc.h"
58 #include "selelem.h"
60 namespace gmx
63 /********************************************************************
64 * SelectionParserSymbol
67 /*! \internal \brief
68 * Private implementation class for SelectionParserSymbol.
70 * \ingroup module_selection
72 class SelectionParserSymbol::Impl
74 public:
75 /*! \brief
76 * Initializes a symbol.
78 * \param[in] type Type for the symbol.
79 * \param[in] name Name for the symbol.
81 * The symbol table is responsible for initializing the \a meth_ and
82 * \a var_ members as appropriate.
84 Impl(SymbolType type, const char* name) : name_(name), type_(type), meth_(nullptr) {}
86 //! Name of the symbol.
87 std::string name_;
88 //! Type of the symbol.
89 SymbolType type_;
90 //! Pointer to the method structure (\ref MethodSymbol).
91 gmx_ana_selmethod_t* meth_;
92 //! Pointer to the variable value (\ref VariableSymbol).
93 SelectionTreeElementPointer var_;
96 SelectionParserSymbol::SelectionParserSymbol(Impl* impl) : impl_(impl) {}
98 SelectionParserSymbol::~SelectionParserSymbol() {}
100 const std::string& SelectionParserSymbol::name() const
102 return impl_->name_;
105 SelectionParserSymbol::SymbolType SelectionParserSymbol::type() const
107 return impl_->type_;
110 gmx_ana_selmethod_t* SelectionParserSymbol::methodValue() const
112 GMX_RELEASE_ASSERT(type() == MethodSymbol,
113 "Attempting to get method handle for a non-method symbol");
114 return impl_->meth_;
117 const gmx::SelectionTreeElementPointer& SelectionParserSymbol::variableValue() const
119 GMX_RELEASE_ASSERT(type() == VariableSymbol,
120 "Attempting to get variable value for a non-variable symbol");
121 return impl_->var_;
124 /********************************************************************
125 * SelectionParserSymbolTable::Impl
128 /*! \internal
129 * \brief
130 * Private implementation class for SelectionParserSymbolTable.
132 * All methods in this class may throw std::bad_alloc if out of memory.
134 * \ingroup module_selection
136 class SelectionParserSymbolTable::Impl
138 public:
139 //! Smart pointer type for managing a SelectionParserSymbol.
140 typedef std::unique_ptr<SelectionParserSymbol> SymbolPointer;
141 //! Container type for the list of symbols.
142 typedef std::map<std::string, SymbolPointer> SymbolMap;
144 /*! \brief
145 * Adds a symbol to the symbol list.
147 * \param[in] symbol Symbol to add.
149 void addSymbol(SymbolPointer symbol);
150 //! Adds the reserved symbols to this symbol table.
151 void addReservedSymbols();
152 //! Adds the position symbols to this symbol table.
153 void addPositionSymbols();
155 //! Symbols in this symbol table.
156 SymbolMap symbols_;
159 void SelectionParserSymbolTable::Impl::addSymbol(SymbolPointer symbol)
161 symbols_.insert(std::make_pair(symbol->name(), std::move(symbol)));
164 void SelectionParserSymbolTable::Impl::addReservedSymbols()
166 const char* const sym_reserved[] = { "group", "to", "not", "and", "or",
167 "xor", "yes", "no", "on", "off" };
169 for (int i = 0; i < asize(sym_reserved); ++i)
171 SymbolPointer sym(new SelectionParserSymbol(new SelectionParserSymbol::Impl(
172 SelectionParserSymbol::ReservedSymbol, sym_reserved[i])));
173 addSymbol(std::move(sym));
177 void SelectionParserSymbolTable::Impl::addPositionSymbols()
179 const char* const* postypes = gmx::PositionCalculationCollection::typeEnumValues;
180 for (int i = 0; postypes[i] != nullptr; ++i)
182 SymbolPointer sym(new SelectionParserSymbol(new SelectionParserSymbol::Impl(
183 SelectionParserSymbol::PositionSymbol, postypes[i])));
184 addSymbol(std::move(sym));
188 /********************************************************************
189 * SelectionParserSymbolIterator
192 /*! \internal \brief
193 * Private implementation class for SelectionParserSymbolIterator.
195 * \ingroup module_selection
197 class SelectionParserSymbolIterator::Impl
199 public:
200 //! Shorthand for the underlying iterator type.
201 typedef SelectionParserSymbolTable::Impl::SymbolMap::const_iterator IteratorType;
203 /*! \brief
204 * Constructs an end iterator.
206 * \param[in] end Iterator to the end of the iterated container.
208 explicit Impl(IteratorType end) : iter_(end), end_(end) {}
209 /*! \brief
210 * Constructs an iterator.
212 * \param[in] iter Iterator to the current symbol.
213 * \param[in] end Iterator to the end of the iterated container.
215 Impl(IteratorType iter, IteratorType end) : iter_(iter), end_(end) {}
217 //! Underlying iterator to the symbol container.
218 IteratorType iter_;
219 //! End of the symbol container being iterated.
220 IteratorType end_;
223 SelectionParserSymbolIterator::SelectionParserSymbolIterator(Impl* impl) : impl_(impl) {}
225 SelectionParserSymbolIterator::SelectionParserSymbolIterator(const SelectionParserSymbolIterator& other) :
226 impl_(new Impl(*other.impl_))
230 SelectionParserSymbolIterator::~SelectionParserSymbolIterator() {}
232 SelectionParserSymbolIterator& SelectionParserSymbolIterator::operator=(const SelectionParserSymbolIterator& other)
234 impl_.reset(new Impl(*other.impl_));
235 return *this;
238 bool SelectionParserSymbolIterator::operator==(const SelectionParserSymbolIterator& other) const
240 return impl_->iter_ == other.impl_->iter_;
243 const SelectionParserSymbol& SelectionParserSymbolIterator::operator*() const
245 return *impl_->iter_->second;
248 SelectionParserSymbolIterator& SelectionParserSymbolIterator::operator++()
250 SelectionParserSymbol::SymbolType type = impl_->iter_->second->type();
253 ++impl_->iter_;
254 } while (impl_->iter_ != impl_->end_ && impl_->iter_->second->type() != type);
255 return *this;
258 /********************************************************************
259 * SelectionParserSymbolTable
262 SelectionParserSymbolTable::SelectionParserSymbolTable() : impl_(new Impl)
264 impl_->addReservedSymbols();
265 impl_->addPositionSymbols();
268 SelectionParserSymbolTable::~SelectionParserSymbolTable() {}
270 const SelectionParserSymbol* SelectionParserSymbolTable::findSymbol(const std::string& name) const
272 Impl::SymbolMap::const_iterator sym = impl_->symbols_.lower_bound(name);
273 if (sym == impl_->symbols_.end())
275 return nullptr;
277 if (sym->second->name() == name)
279 return sym->second.get();
281 return nullptr;
284 SelectionParserSymbolIterator SelectionParserSymbolTable::beginIterator(SelectionParserSymbol::SymbolType type) const
286 Impl::SymbolMap::const_iterator sym;
287 Impl::SymbolMap::const_iterator end = impl_->symbols_.end();
288 for (sym = impl_->symbols_.begin(); sym != end; ++sym)
290 if (sym->second->type() == type)
292 return SelectionParserSymbolIterator(new SelectionParserSymbolIterator::Impl(sym, end));
295 return endIterator();
298 SelectionParserSymbolIterator SelectionParserSymbolTable::endIterator() const
300 return SelectionParserSymbolIterator(new SelectionParserSymbolIterator::Impl(impl_->symbols_.end()));
303 void SelectionParserSymbolTable::addVariable(const char* name, const gmx::SelectionTreeElementPointer& sel)
305 // In the current parser implementation, a syntax error is produced before
306 // this point is reached, but the check is here for robustness.
307 Impl::SymbolMap::const_iterator other = impl_->symbols_.find(name);
308 if (other != impl_->symbols_.end())
310 if (other->second->type() == SelectionParserSymbol::VariableSymbol)
312 GMX_THROW(InvalidInputError(formatString("Reassigning variable '%s' is not supported", name)));
314 else
316 GMX_THROW(InvalidInputError(
317 formatString("Variable name '%s' conflicts with a reserved keyword", name)));
320 Impl::SymbolPointer sym(new SelectionParserSymbol(
321 new SelectionParserSymbol::Impl(SelectionParserSymbol::VariableSymbol, name)));
322 sym->impl_->var_ = sel;
323 impl_->addSymbol(std::move(sym));
326 void SelectionParserSymbolTable::addMethod(const char* name, gmx_ana_selmethod_t* method)
328 if (impl_->symbols_.find(name) != impl_->symbols_.end())
330 GMX_THROW(APIError(formatString("Method name '%s' conflicts with another symbol", name)));
332 Impl::SymbolPointer sym(new SelectionParserSymbol(
333 new SelectionParserSymbol::Impl(SelectionParserSymbol::MethodSymbol, name)));
334 sym->impl_->meth_ = method;
335 impl_->addSymbol(std::move(sym));
338 } // namespace gmx