Reduce preprocessor dependency in resourcedivision.cpp
[gromacs.git] / src / gromacs / selection / symrec.cpp
blob1cbc5df4794de7171d0f4f8f3452782735c5ecdc
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2009,2010,2011,2012,2013,2014,2015,2017, 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 /*! \internal \file
36 * \brief
37 * Implements classes in symrec.h.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_selection
42 #include "gmxpre.h"
44 #include "symrec.h"
46 #include <map>
47 #include <memory>
48 #include <string>
49 #include <utility>
51 #include "gromacs/utility/arraysize.h"
52 #include "gromacs/utility/exceptions.h"
53 #include "gromacs/utility/gmxassert.h"
54 #include "gromacs/utility/stringutil.h"
56 #include "poscalc.h"
57 #include "selelem.h"
59 namespace gmx
62 /********************************************************************
63 * SelectionParserSymbol
66 /*! \internal \brief
67 * Private implementation class for SelectionParserSymbol.
69 * \ingroup module_selection
71 class SelectionParserSymbol::Impl
73 public:
74 /*! \brief
75 * Initializes a symbol.
77 * \param[in] type Type for the symbol.
78 * \param[in] name Name for the symbol.
80 * The symbol table is responsible for initializing the \a meth_ and
81 * \a var_ members as appropriate.
83 Impl(SymbolType type, const char *name)
84 : name_(name), type_(type), meth_(nullptr)
88 //! Name of the symbol.
89 std::string name_;
90 //! Type of the symbol.
91 SymbolType type_;
92 //! Pointer to the method structure (\ref MethodSymbol).
93 gmx_ana_selmethod_t *meth_;
94 //! Pointer to the variable value (\ref VariableSymbol).
95 SelectionTreeElementPointer var_;
98 SelectionParserSymbol::SelectionParserSymbol(Impl *impl)
99 : impl_(impl)
103 SelectionParserSymbol::~SelectionParserSymbol()
107 const std::string &
108 SelectionParserSymbol::name() const
110 return impl_->name_;
113 SelectionParserSymbol::SymbolType
114 SelectionParserSymbol::type() const
116 return impl_->type_;
119 gmx_ana_selmethod_t *
120 SelectionParserSymbol::methodValue() const
122 GMX_RELEASE_ASSERT(type() == MethodSymbol,
123 "Attempting to get method handle for a non-method symbol");
124 return impl_->meth_;
127 const gmx::SelectionTreeElementPointer &
128 SelectionParserSymbol::variableValue() const
130 GMX_RELEASE_ASSERT(type() == VariableSymbol,
131 "Attempting to get variable value for a non-variable symbol");
132 return impl_->var_;
135 /********************************************************************
136 * SelectionParserSymbolTable::Impl
139 /*! \internal
140 * \brief
141 * Private implementation class for SelectionParserSymbolTable.
143 * All methods in this class may throw std::bad_alloc if out of memory.
145 * \ingroup module_selection
147 class SelectionParserSymbolTable::Impl
149 public:
150 //! Smart pointer type for managing a SelectionParserSymbol.
151 typedef std::unique_ptr<SelectionParserSymbol>
152 SymbolPointer;
153 //! Container type for the list of symbols.
154 typedef std::map<std::string, SymbolPointer> SymbolMap;
156 /*! \brief
157 * Adds a symbol to the symbol list.
159 * \param[in] symbol Symbol to add.
161 void addSymbol(SymbolPointer symbol);
162 //! Adds the reserved symbols to this symbol table.
163 void addReservedSymbols();
164 //! Adds the position symbols to this symbol table.
165 void addPositionSymbols();
167 //! Symbols in this symbol table.
168 SymbolMap symbols_;
171 void
172 SelectionParserSymbolTable::Impl::addSymbol(SymbolPointer symbol)
174 symbols_.insert(std::make_pair(symbol->name(), std::move(symbol)));
177 void
178 SelectionParserSymbolTable::Impl::addReservedSymbols()
180 const char *const sym_reserved[] = {
181 "group",
182 "to",
183 "not",
184 "and",
185 "or",
186 "xor",
187 "yes",
188 "no",
189 "on",
190 "off"
193 for (size_t i = 0; i < asize(sym_reserved); ++i)
195 SymbolPointer sym(new SelectionParserSymbol(
196 new SelectionParserSymbol::Impl(
197 SelectionParserSymbol::ReservedSymbol, sym_reserved[i])));
198 addSymbol(std::move(sym));
202 void
203 SelectionParserSymbolTable::Impl::addPositionSymbols()
205 const char *const *postypes
206 = gmx::PositionCalculationCollection::typeEnumValues;
207 for (int i = 0; postypes[i] != nullptr; ++i)
209 SymbolPointer sym(new SelectionParserSymbol(
210 new SelectionParserSymbol::Impl(
211 SelectionParserSymbol::PositionSymbol, postypes[i])));
212 addSymbol(std::move(sym));
216 /********************************************************************
217 * SelectionParserSymbolIterator
220 /*! \internal \brief
221 * Private implementation class for SelectionParserSymbolIterator.
223 * \ingroup module_selection
225 class SelectionParserSymbolIterator::Impl
227 public:
228 //! Shorthand for the underlying iterator type.
229 typedef SelectionParserSymbolTable::Impl::SymbolMap::const_iterator
230 IteratorType;
232 /*! \brief
233 * Constructs an end iterator.
235 * \param[in] end Iterator to the end of the iterated container.
237 explicit Impl(IteratorType end)
238 : iter_(end), end_(end)
241 /*! \brief
242 * Constructs an iterator.
244 * \param[in] iter Iterator to the current symbol.
245 * \param[in] end Iterator to the end of the iterated container.
247 Impl(IteratorType iter, IteratorType end)
248 : iter_(iter), end_(end)
252 //! Underlying iterator to the symbol container.
253 IteratorType iter_;
254 //! End of the symbol container being iterated.
255 IteratorType end_;
258 SelectionParserSymbolIterator::SelectionParserSymbolIterator(Impl *impl)
259 : impl_(impl)
263 SelectionParserSymbolIterator::SelectionParserSymbolIterator(
264 const SelectionParserSymbolIterator &other)
265 : impl_(new Impl(*other.impl_))
269 SelectionParserSymbolIterator::~SelectionParserSymbolIterator()
273 SelectionParserSymbolIterator &SelectionParserSymbolIterator::operator=(
274 const SelectionParserSymbolIterator &other)
276 impl_.reset(new Impl(*other.impl_));
277 return *this;
280 bool SelectionParserSymbolIterator::operator==(
281 const SelectionParserSymbolIterator &other) const
283 return impl_->iter_ == other.impl_->iter_;
286 const SelectionParserSymbol &SelectionParserSymbolIterator::operator*() const
288 return *impl_->iter_->second;
291 SelectionParserSymbolIterator &SelectionParserSymbolIterator::operator++()
293 SelectionParserSymbol::SymbolType type = impl_->iter_->second->type();
296 ++impl_->iter_;
298 while (impl_->iter_ != impl_->end_ && impl_->iter_->second->type() != type);
299 return *this;
302 /********************************************************************
303 * SelectionParserSymbolTable
306 SelectionParserSymbolTable::SelectionParserSymbolTable()
307 : impl_(new Impl)
309 impl_->addReservedSymbols();
310 impl_->addPositionSymbols();
313 SelectionParserSymbolTable::~SelectionParserSymbolTable()
317 const SelectionParserSymbol *
318 SelectionParserSymbolTable::findSymbol(const std::string &name) const
320 Impl::SymbolMap::const_iterator sym = impl_->symbols_.lower_bound(name);
321 if (sym == impl_->symbols_.end())
323 return nullptr;
325 if (sym->second->name() == name)
327 return sym->second.get();
329 return nullptr;
332 SelectionParserSymbolIterator
333 SelectionParserSymbolTable::beginIterator(SelectionParserSymbol::SymbolType type) const
335 Impl::SymbolMap::const_iterator sym;
336 Impl::SymbolMap::const_iterator end = impl_->symbols_.end();
337 for (sym = impl_->symbols_.begin(); sym != end; ++sym)
339 if (sym->second->type() == type)
341 return SelectionParserSymbolIterator(
342 new SelectionParserSymbolIterator::Impl(sym, end));
345 return endIterator();
348 SelectionParserSymbolIterator
349 SelectionParserSymbolTable::endIterator() const
351 return SelectionParserSymbolIterator(
352 new SelectionParserSymbolIterator::Impl(impl_->symbols_.end()));
355 void
356 SelectionParserSymbolTable::addVariable(const char *name,
357 const gmx::SelectionTreeElementPointer &sel)
359 // In the current parser implementation, a syntax error is produced before
360 // this point is reached, but the check is here for robustness.
361 Impl::SymbolMap::const_iterator other = impl_->symbols_.find(name);
362 if (other != impl_->symbols_.end())
364 if (other->second->type() == SelectionParserSymbol::VariableSymbol)
366 GMX_THROW(InvalidInputError(
367 formatString("Reassigning variable '%s' is not supported",
368 name)));
370 else
372 GMX_THROW(InvalidInputError(
373 formatString("Variable name '%s' conflicts with a reserved keyword",
374 name)));
377 Impl::SymbolPointer sym(new SelectionParserSymbol(
378 new SelectionParserSymbol::Impl(
379 SelectionParserSymbol::VariableSymbol, name)));
380 sym->impl_->var_ = sel;
381 impl_->addSymbol(std::move(sym));
384 void
385 SelectionParserSymbolTable::addMethod(const char *name,
386 gmx_ana_selmethod_t *method)
388 if (impl_->symbols_.find(name) != impl_->symbols_.end())
390 GMX_THROW(APIError(
391 formatString("Method name '%s' conflicts with another symbol",
392 name)));
394 Impl::SymbolPointer sym(new SelectionParserSymbol(
395 new SelectionParserSymbol::Impl(
396 SelectionParserSymbol::MethodSymbol, name)));
397 sym->impl_->meth_ = method;
398 impl_->addSymbol(std::move(sym));
401 } // namespace gmx