Split off nbnxn GPU timing and staging reduction
[gromacs.git] / src / gromacs / options / basicoptions.cpp
blob57e691f74f0dfb08aebb7eed27a937d5d569c241
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2010,2011,2012,2013,2014,2015,2016,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 basicoptions.h and basicoptionstorage.h.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_options
42 #include "gmxpre.h"
44 #include "basicoptions.h"
46 #include <cerrno>
47 #include <cstdio>
48 #include <cstdlib>
50 #include <algorithm>
51 #include <limits>
52 #include <string>
53 #include <vector>
55 #include "gromacs/utility/cstringutil.h"
56 #include "gromacs/utility/exceptions.h"
57 #include "gromacs/utility/strconvert.h"
58 #include "gromacs/utility/stringutil.h"
60 #include "basicoptionstorage.h"
62 namespace
65 /*! \brief
66 * Expands a single value to a vector by copying the value.
68 * \tparam ValueType Type of values to process.
69 * \param[in] length Length of the resulting vector.
70 * \param[in,out] values Values to process.
71 * \throws std::bad_alloc if out of memory.
72 * \throws InvalidInputError if \p values has an invalid number of values.
74 * \p values should have 0, 1, or \p length values.
75 * If \p values has 1 value, it is expanded such that it has \p length
76 * identical values. In other valid cases, nothing is done.
78 * \ingroup module_options
80 template <typename ValueType>
81 void expandVector(size_t length, std::vector<ValueType> *values)
83 if (length > 0 && !values->empty() && values->size() != length)
85 if (values->size() != 1)
87 GMX_THROW(gmx::InvalidInputError(gmx::formatString(
88 "Expected 1 or %d values, got %d", length, values->size())));
90 const ValueType &value = (*values)[0];
91 values->resize(length, value);
95 /*! \brief
96 * Finds an enumerated value from the list of allowed values.
98 * \param[in] allowedValues List of allowed values.
99 * \param[in] value Value to search for.
100 * \throws gmx::InvalidInputError if \p value does not match anything in
101 * \p allowedValues.
102 * \returns Iterator to the found value.
104 * \ingroup module_options
106 std::vector<std::string>::const_iterator
107 findEnumValue(const std::vector<std::string> &allowedValues,
108 const std::string &value)
110 std::vector<std::string>::const_iterator i;
111 std::vector<std::string>::const_iterator match = allowedValues.end();
112 for (i = allowedValues.begin(); i != allowedValues.end(); ++i)
114 // TODO: Case independence.
115 if (gmx::startsWith(*i, value))
117 if (match == allowedValues.end() || i->size() < match->size())
119 match = i;
123 if (match == allowedValues.end())
125 GMX_THROW(gmx::InvalidInputError("Invalid value: " + value));
127 return match;
130 } // namespace
132 namespace gmx
135 /********************************************************************
136 * BooleanOptionStorage
139 std::string BooleanOptionStorage::formatSingleValue(const bool &value) const
141 return value ? "yes" : "no";
144 void BooleanOptionStorage::initConverter(ConverterType *converter)
146 converter->addConverter<std::string>(&fromStdString<bool>);
149 /********************************************************************
150 * BooleanOptionInfo
153 BooleanOptionInfo::BooleanOptionInfo(BooleanOptionStorage *option)
154 : OptionInfo(option)
158 const BooleanOptionStorage &BooleanOptionInfo::option() const
160 return static_cast<const BooleanOptionStorage &>(OptionInfo::option());
163 bool BooleanOptionInfo::defaultValue() const
165 return option().defaultValue();
168 /********************************************************************
169 * BooleanOption
172 AbstractOptionStorage *
173 BooleanOption::createStorage(const OptionManagerContainer & /*managers*/) const
175 return new BooleanOptionStorage(*this);
179 /********************************************************************
180 * IntegerOptionStorage
183 std::string IntegerOptionStorage::formatSingleValue(const int &value) const
185 return toString(value);
188 void IntegerOptionStorage::initConverter(ConverterType *converter)
190 converter->addConverter<std::string>(&fromStdString<int>);
193 void IntegerOptionStorage::processSetValues(ValueList *values)
195 if (isVector())
197 expandVector(maxValueCount(), values);
201 /********************************************************************
202 * IntegerOptionInfo
205 IntegerOptionInfo::IntegerOptionInfo(IntegerOptionStorage *option)
206 : OptionInfo(option)
210 /********************************************************************
211 * IntegerOption
214 AbstractOptionStorage *
215 IntegerOption::createStorage(const OptionManagerContainer & /*managers*/) const
217 return new IntegerOptionStorage(*this);
221 /********************************************************************
222 * Int64OptionStorage
225 std::string Int64OptionStorage::formatSingleValue(const gmx_int64_t &value) const
227 return toString(value);
230 void Int64OptionStorage::initConverter(ConverterType *converter)
232 converter->addConverter<std::string>(&fromStdString<gmx_int64_t>);
235 /********************************************************************
236 * Int64OptionInfo
239 Int64OptionInfo::Int64OptionInfo(Int64OptionStorage *option)
240 : OptionInfo(option)
244 /********************************************************************
245 * Int64Option
248 AbstractOptionStorage *
249 Int64Option::createStorage(const OptionManagerContainer & /*managers*/) const
251 return new Int64OptionStorage(*this);
255 /********************************************************************
256 * DoubleOptionStorage
259 DoubleOptionStorage::DoubleOptionStorage(const DoubleOption &settings)
260 : MyBase(settings), info_(this), bTime_(settings.bTime_), factor_(1.0)
264 std::string DoubleOptionStorage::typeString() const
266 return isVector() ? "vector" : (isTime() ? "time" : "real");
269 std::string DoubleOptionStorage::formatSingleValue(const double &value) const
271 return toString(value / factor_);
274 void DoubleOptionStorage::initConverter(ConverterType *converter)
276 converter->addConverter<std::string>(&fromStdString<double>);
277 converter->addCastConversion<float>();
280 double DoubleOptionStorage::processValue(const double &value) const
282 // TODO: Consider testing for overflow when scaling with factor_.
283 return value * factor_;
286 void DoubleOptionStorage::processSetValues(ValueList *values)
288 if (isVector())
290 expandVector(maxValueCount(), values);
294 void DoubleOptionStorage::setScaleFactor(double factor)
296 GMX_RELEASE_ASSERT(factor > 0.0, "Invalid scaling factor");
297 if (!hasFlag(efOption_HasDefaultValue))
299 double scale = factor / factor_;
300 for (double &value : values())
302 value *= scale;
305 factor_ = factor;
308 /********************************************************************
309 * DoubleOptionInfo
312 DoubleOptionInfo::DoubleOptionInfo(DoubleOptionStorage *option)
313 : OptionInfo(option)
317 DoubleOptionStorage &DoubleOptionInfo::option()
319 return static_cast<DoubleOptionStorage &>(OptionInfo::option());
322 const DoubleOptionStorage &DoubleOptionInfo::option() const
324 return static_cast<const DoubleOptionStorage &>(OptionInfo::option());
327 bool DoubleOptionInfo::isTime() const
329 return option().isTime();
332 void DoubleOptionInfo::setScaleFactor(double factor)
334 option().setScaleFactor(factor);
337 /********************************************************************
338 * DoubleOption
341 AbstractOptionStorage *
342 DoubleOption::createStorage(const OptionManagerContainer & /*managers*/) const
344 return new DoubleOptionStorage(*this);
348 /********************************************************************
349 * FloatOptionStorage
352 FloatOptionStorage::FloatOptionStorage(const FloatOption &settings)
353 : MyBase(settings), info_(this), bTime_(settings.bTime_), factor_(1.0)
357 std::string FloatOptionStorage::typeString() const
359 return isVector() ? "vector" : (isTime() ? "time" : "real");
362 std::string FloatOptionStorage::formatSingleValue(const float &value) const
364 return toString(value / factor_);
367 void FloatOptionStorage::initConverter(ConverterType *converter)
369 converter->addConverter<std::string>(&fromStdString<float>);
370 converter->addCastConversion<double>();
373 float FloatOptionStorage::processValue(const float &value) const
375 // TODO: Consider testing for overflow when scaling with factor_.
376 return value * factor_;
379 void FloatOptionStorage::processSetValues(ValueList *values)
381 if (isVector())
383 expandVector(maxValueCount(), values);
387 void FloatOptionStorage::setScaleFactor(double factor)
389 GMX_RELEASE_ASSERT(factor > 0.0, "Invalid scaling factor");
390 if (!hasFlag(efOption_HasDefaultValue))
392 float scale = factor / factor_;
393 for (float &value : values())
395 value *= scale;
398 factor_ = factor;
401 /********************************************************************
402 * FloatOptionInfo
405 FloatOptionInfo::FloatOptionInfo(FloatOptionStorage *option)
406 : OptionInfo(option)
410 FloatOptionStorage &FloatOptionInfo::option()
412 return static_cast<FloatOptionStorage &>(OptionInfo::option());
415 const FloatOptionStorage &FloatOptionInfo::option() const
417 return static_cast<const FloatOptionStorage &>(OptionInfo::option());
420 bool FloatOptionInfo::isTime() const
422 return option().isTime();
425 void FloatOptionInfo::setScaleFactor(double factor)
427 option().setScaleFactor(factor);
430 /********************************************************************
431 * FloatOption
434 AbstractOptionStorage *
435 FloatOption::createStorage(const OptionManagerContainer & /*managers*/) const
437 return new FloatOptionStorage(*this);
441 /********************************************************************
442 * StringOptionStorage
445 StringOptionStorage::StringOptionStorage(const StringOption &settings)
446 : MyBase(settings), info_(this)
448 if (settings.defaultEnumIndex_ >= 0 && settings.enumValues_ == nullptr)
450 GMX_THROW(APIError("Cannot set default enum index without enum values"));
452 if (settings.enumValues_ != nullptr)
454 int count = settings.enumValuesCount_;
455 if (count < 0)
457 count = 0;
458 while (settings.enumValues_[count] != nullptr)
460 ++count;
463 for (int i = 0; i < count; ++i)
465 if (settings.enumValues_[i] == nullptr)
467 GMX_THROW(APIError("Enumeration value cannot be NULL"));
469 allowed_.emplace_back(settings.enumValues_[i]);
471 if (settings.defaultEnumIndex_ >= 0)
473 if (settings.defaultEnumIndex_ >= count)
475 GMX_THROW(APIError("Default enumeration index is out of range"));
477 const std::string *defaultValue = settings.defaultValue();
478 if (defaultValue != nullptr && *defaultValue != allowed_[settings.defaultEnumIndex_])
480 GMX_THROW(APIError("Conflicting default values"));
482 setDefaultValue(allowed_[settings.defaultEnumIndex_]);
487 std::string StringOptionStorage::formatExtraDescription() const
489 std::string result;
490 if (!allowed_.empty())
492 result.append(": ");
493 result.append(joinStrings(allowed_, ", "));
495 return result;
498 std::string StringOptionStorage::formatSingleValue(const std::string &value) const
500 return value;
503 void StringOptionStorage::initConverter(ConverterType * /*converter*/)
507 std::string StringOptionStorage::processValue(const std::string &value) const
509 if (allowed_.size() > 0)
511 return *findEnumValue(this->allowed_, value);
513 return value;
516 /********************************************************************
517 * StringOptionInfo
520 StringOptionInfo::StringOptionInfo(StringOptionStorage *option)
521 : OptionInfo(option)
525 const StringOptionStorage &StringOptionInfo::option() const
527 return static_cast<const StringOptionStorage &>(OptionInfo::option());
530 bool StringOptionInfo::isEnumerated() const
532 return !allowedValues().empty();
535 const std::vector<std::string> &StringOptionInfo::allowedValues() const
537 return option().allowedValues();
540 /********************************************************************
541 * StringOption
544 AbstractOptionStorage *
545 StringOption::createStorage(const OptionManagerContainer & /*managers*/) const
547 return new StringOptionStorage(*this);
551 /********************************************************************
552 * EnumOptionStorage
555 EnumOptionStorage::EnumOptionStorage(const AbstractOption &settings,
556 const char *const *enumValues, int count,
557 int defaultValue, int defaultValueIfSet,
558 StorePointer store)
559 : MyBase(settings, std::move(store)), info_(this)
561 if (enumValues == nullptr)
563 GMX_THROW(APIError("Allowed values must be provided to EnumOption"));
566 if (count < 0)
568 count = 0;
569 while (enumValues[count] != nullptr)
571 ++count;
574 for (int i = 0; i < count; ++i)
576 if (enumValues[i] == nullptr)
578 GMX_THROW(APIError("Enumeration value cannot be NULL"));
580 allowed_.emplace_back(enumValues[i]);
583 GMX_ASSERT(defaultValue < count, "Default enumeration value is out of range");
584 GMX_ASSERT(defaultValueIfSet < count, "Default enumeration value is out of range");
585 setFlag(efOption_HasDefaultValue);
586 if (defaultValue >= 0)
588 setDefaultValue(defaultValue);
590 if (defaultValueIfSet >= 0)
592 setDefaultValueIfSet(defaultValueIfSet);
596 std::string EnumOptionStorage::formatExtraDescription() const
598 std::string result;
599 result.append(": ");
600 result.append(joinStrings(allowed_, ", "));
601 return result;
604 std::string EnumOptionStorage::formatSingleValue(const int &value) const
606 if (value < 0 || value >= static_cast<int>(allowed_.size()))
608 return std::string();
610 return allowed_[value];
613 Variant EnumOptionStorage::normalizeValue(const int &value) const
615 return Variant::create<std::string>(formatSingleValue(value));
618 void EnumOptionStorage::initConverter(ConverterType *converter)
620 converter->addConverter<std::string>(
621 [this] (const std::string &value)
623 return findEnumValue(this->allowed_, value) - this->allowed_.begin();
627 /********************************************************************
628 * EnumOptionInfo
631 EnumOptionInfo::EnumOptionInfo(EnumOptionStorage *option)
632 : OptionInfo(option)
636 const EnumOptionStorage &EnumOptionInfo::option() const
638 return static_cast<const EnumOptionStorage &>(OptionInfo::option());
641 const std::vector<std::string> &EnumOptionInfo::allowedValues() const
643 return option().allowedValues();
646 /********************************************************************
647 * EnumOption helpers
650 namespace internal
653 //! \cond internal
654 AbstractOptionStorage *
655 createEnumOptionStorage(const AbstractOption &option,
656 const char *const *enumValues, int count,
657 int defaultValue, int defaultValueIfSet,
658 IOptionValueStore<int> *store)
660 std::unique_ptr<IOptionValueStore<int> > storePtr(store);
661 return new EnumOptionStorage(option, enumValues, count, defaultValue,
662 defaultValueIfSet, move(storePtr));
664 //! \endcond
666 } // namespace internal
668 } // namespace gmx