Fix clang 10 -Wimplicit-int-float-conversion warnings
[gromacs.git] / src / gromacs / random / uniformrealdistribution.h
blobbc7cb61405354569caf0207b387153f16db3722f
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2015,2016,2018,2019,2020, 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.
36 /*! \file
37 * \brief The uniform real distribution
39 * Portable version of the uniform real that generates the same sequence
40 * on all platforms. Since stdlibc++ and libc++ provide different sequences
41 * we prefer this one so unit tests produce the same values on all platforms.
43 * \author Erik Lindahl <erik.lindahl@gmail.com>
44 * \inpublicapi
45 * \ingroup module_random
48 #ifndef GMX_RANDOM_UNIFORMREALDISTRIBUTION_H
49 #define GMX_RANDOM_UNIFORMREALDISTRIBUTION_H
51 #include <cmath>
53 #include <algorithm>
54 #include <limits>
55 #include <type_traits>
57 #include "gromacs/math/functions.h"
58 #include "gromacs/utility/basedefinitions.h"
59 #include "gromacs/utility/classhelpers.h"
60 #include "gromacs/utility/gmxassert.h"
61 #include "gromacs/utility/real.h"
64 * The portable version of the uniform real distribution (to make sure we get
65 * the same values on all platforms) has been modified from the LLVM libcxx
66 * headers, distributed under the MIT license:
68 * Copyright (c) The LLVM compiler infrastructure
70 * Permission is hereby granted, free of charge, to any person obtaining a copy
71 * of this software and associated documentation files (the "Software"), to deal
72 * in the Software without restriction, including without limitation the rights
73 * to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
74 * copies of the Software, and to permit persons to whom the Software is
75 * furnished to do so, subject to the following conditions:
77 * The above copyright notice and this permission notice shall be included in
78 * all copies or substantial portions of the Software.
80 * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
81 * IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
82 * FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
83 * AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
84 * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
85 * OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN
86 * THE SOFTWARE.
89 namespace gmx
92 /*! \brief Generate a floating-point value with specified number of random bits
94 * \tparam RealType Floating-point type to generate
95 * \tparam Bits Number of random bits to generate
96 * \tparam Rng Random number generator class
98 * \param g Random number generator to use
100 * This implementation avoids the bug in libc++ and stdlibc++ (which is due
101 * to the C++ standard being unclear) where 1.0 can be returned occasionally.
104 template<class RealType = real, unsigned int Bits, class Rng>
105 RealType generateCanonical(Rng& g)
107 // No point in using more bits than fit in RealType
108 const uint64_t digits = std::numeric_limits<RealType>::digits;
109 const uint64_t realBits = std::min(digits, static_cast<uint64_t>(Bits));
110 const uint64_t range = Rng::max() - Rng::min() + uint64_t(1);
111 uint64_t log2R = (range == 0) ? std::numeric_limits<uint64_t>::digits : log2I(range);
112 uint64_t k = realBits / log2R + (realBits % log2R != 0) + (realBits == 0);
113 // Note that Rng::max and Rng::min are typically an integer type.
114 // Only unsigned integer types can express the range using the
115 // same type. Converting to RealType before computing the range
116 // would work but we have no need for that.
117 static_assert(std::is_unsigned<decltype(Rng::max())>::value
118 && std::is_unsigned<decltype(Rng::min())>::value,
119 "Rng::max and Rng::min must be unsigned");
120 RealType r = RealType(Rng::max() - Rng::min()) + RealType(1);
121 RealType s = g() - Rng::min();
122 RealType base = r;
123 RealType result;
125 for (uint64_t i = 1; i < k; ++i)
127 s += RealType(g() - Rng::min()) * base;
128 base *= r;
130 result = s / base;
132 // This implementation is specified by the C++ standard, but unfortunately it
133 // has a bug where 1.0 can be generated occasionally due to the limited
134 // precision of floating point, while 0.0 is only generated half as often as
135 // it should. We "solve" both these issues by swapping 1.0 for 0.0 when it happens.
137 // See:
138 // https://llvm.org/bugs/show_bug.cgi?id=18767
139 // https://gcc.gnu.org/bugzilla/show_bug.cgi?id=63176
141 // Note that we prefer not to use the gcc 'fix' of looping until the result
142 // is smaller than 1.0, since that breaks the strict specification of the
143 // number of times the rng will be called.
145 // This can only happen when we ask for the same number of bits that fit
146 // in RealType, so by checking for that we avoid the extra code in all other
147 // cases. If you are worried about it: Use RealType=double with 32 bits.
149 if (realBits == digits && result == 1.0)
151 result = 0.0;
153 return result;
157 /*! \brief Uniform real distribution
159 * The C++ standard library does provide this distribution, but even
160 * though they all sample from the correct distribution different standard
161 * library implementations appear to return different sequences of numbers
162 * for the same random number generator. To make it easier to use GROMACS
163 * unit tests that depend on random numbers we have our own implementation.
165 * \tparam RealType Floating-point type, real by default in GROMACS.
167 template<class RealType = real>
168 class UniformRealDistribution
170 public:
171 /*! \brief Type of values returned */
172 typedef RealType result_type;
174 /*! \brief Uniform real distribution parameters */
175 class param_type
177 /*! \brief Lower end of range (inclusive) */
178 result_type a_;
179 /*! \brief Upper end of range (exclusive) */
180 result_type b_;
182 public:
183 /*! \brief Reference back to the distribution class */
184 typedef UniformRealDistribution distribution_type;
186 /*! \brief Construct parameter block
188 * \param a Lower end of range (inclusive)
189 * \param b Upper end of range (exclusive)
191 explicit param_type(result_type a = 0.0, result_type b = 1.0) : a_(a), b_(b)
193 GMX_RELEASE_ASSERT(a < b, "The uniform real distribution requires a<b");
196 /*! \brief Return first parameter */
197 result_type a() const { return a_; }
198 /*! \brief Return second parameter */
199 result_type b() const { return b_; }
201 /*! \brief True if two parameter sets will return the same uniform real distribution.
203 * \param x Instance to compare with.
205 bool operator==(const param_type& x) const { return a_ == x.a_ && b_ == x.b_; }
207 /*! \brief True if two parameter sets will return different uniform real distributions
209 * \param x Instance to compare with.
211 bool operator!=(const param_type& x) const { return !operator==(x); }
214 public:
215 /*! \brief Construct new distribution with given floating-point parameters.
217 * \param a Lower end of range (inclusive)
218 * \param b Upper end of range (exclusive)
220 explicit UniformRealDistribution(result_type a = 0.0, result_type b = 1.0) :
221 param_(param_type(a, b))
225 /*! \brief Construct new distribution from parameter class
227 * \param param Parameter class as defined inside gmx::UniformRealDistribution.
229 explicit UniformRealDistribution(const param_type& param) : param_(param) {}
231 /*! \brief Flush all internal saved values */
232 void reset() {}
234 /*! \brief Return values from uniform real distribution with internal parameters
236 * \tparam Rng Random engine class
238 * \param g Random engine
240 template<class Rng>
241 result_type operator()(Rng& g)
243 return (*this)(g, param_);
246 /*! \brief Return value from uniform real distribution with given parameters
248 * \tparam Rng Random engine class
250 * \param g Random engine
251 * \param param Parameters to use
253 template<class Rng>
254 result_type operator()(Rng& g, const param_type& param)
256 result_type r = generateCanonical<RealType, std::numeric_limits<RealType>::digits>(g);
257 return (param.b() - param.a()) * r + param.a();
260 /*! \brief Return the lower range uniform real distribution */
261 result_type a() const { return param_.a(); }
263 /*! \brief Return the upper range of the uniform real distribution */
264 result_type b() const { return param_.b(); }
266 /*! \brief Return the full parameter class of the uniform real distribution */
267 param_type param() const { return param_; }
269 /*! \brief Smallest value that can be returned from uniform real distribution */
270 result_type min() const { return a(); }
272 /*! \brief Largest value that can be returned from uniform real distribution */
273 result_type max() const { return b(); }
275 /*! \brief True if two uniform real distributions will produce the same values.
277 * \param x Instance to compare with.
279 bool operator==(const UniformRealDistribution& x) const { return param_ == x.param_; }
281 /*! \brief True if two uniform real distributions will produce different values.
283 * \param x Instance to compare with.
285 bool operator!=(const UniformRealDistribution& x) const { return !operator==(x); }
287 private:
288 /*! \brief Internal value for parameters, can be overridden at generation time. */
289 param_type param_;
291 GMX_DISALLOW_COPY_AND_ASSIGN(UniformRealDistribution);
294 } // namespace gmx
296 #endif // GMX_RANDOM_UNIFORMREALDISTRIBUTION_H