Write TPR body as opaque XDR data in big-endian format
[gromacs.git] / src / gromacs / utility / bitmask.h
blobeddec13736eee67b1ff75ff45ddb445e5d63e339
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2018,2019, 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 /*! \libinternal \file
36 * \brief
37 * Declares gmx_bitmask_t and associated functions
39 * \author Roland Schulz <roland@utk.edu>
40 * \inlibraryapi
41 * \ingroup module_utility
44 #ifndef GMX_MDLIB_BITMASK_H
45 #define GMX_MDLIB_BITMASK_H
47 #include "config.h" /* for GMX_MAX_OPENMP_THREADS */
49 #include <string.h>
51 #include <algorithm>
52 #include <array>
54 #include "gromacs/utility/basedefinitions.h"
55 #include "gromacs/utility/stringutil.h"
57 /*! \brief Size of bitmask. Has to be 32 or multiple of 64. */
58 #ifndef BITMASK_SIZE
59 # define BITMASK_SIZE GMX_OPENMP_MAX_THREADS
60 #endif
62 #if BITMASK_SIZE != 32 && BITMASK_SIZE % 64 != 0
63 # error BITMASK_SIZE has to be 32 or a multiple of 64.
64 #endif
66 #if BITMASK_SIZE <= 64 || defined DOXYGEN
67 # if BITMASK_SIZE == 32
68 typedef uint32_t gmx_bitmask_t;
69 # else
70 typedef uint64_t gmx_bitmask_t; /**< bitmask type */
71 # endif
73 /*! \brief Initialize all bits to 0 */
74 inline static void bitmask_clear(gmx_bitmask_t* m)
76 *m = 0;
79 /*! \brief Set bit at position b to 1. */
80 inline static void bitmask_set_bit(gmx_bitmask_t* m, int b)
82 *m |= gmx_bitmask_t(1) << b;
85 /*! \brief Initialize all bits: bit b to 1, others to 0 */
86 inline static void bitmask_init_bit(gmx_bitmask_t* m, int b)
88 *m = gmx_bitmask_t(1) << b;
91 /*! \brief Initialize all bits: all bits below b to 1, others to 0 */
92 inline static void bitmask_init_low_bits(gmx_bitmask_t* m, int b)
94 *m = (gmx_bitmask_t(1) << b) - 1;
97 /*! \brief Test if bit b is set */
98 inline static bool bitmask_is_set(gmx_bitmask_t m, int b)
100 return (m & (gmx_bitmask_t(1) << b)) != 0;
103 /*! \brief Test if both bitmasks have no common bits enabled */
104 inline static bool bitmask_is_disjoint(gmx_bitmask_t a, gmx_bitmask_t b)
106 return (a & b) == 0U;
109 /*! \brief Test if both bitmasks are equal */
110 inline static bool bitmask_is_equal(gmx_bitmask_t a, gmx_bitmask_t b)
112 return a == b;
115 /*! \brief Test if bitmask has no enabled bits */
116 inline static bool bitmask_is_zero(gmx_bitmask_t m)
118 return m == 0U;
121 /*! \brief Set all bits enabled in either mask and write into a */
122 inline static void bitmask_union(gmx_bitmask_t* a, gmx_bitmask_t b)
124 *a |= b;
126 #else
127 # define BITMASK_ALEN (BITMASK_SIZE / 64)
128 using gmx_bitmask_t = std::array<uint64_t, BITMASK_ALEN>;
130 inline static void bitmask_clear(gmx_bitmask_t* m)
132 m->fill(0);
135 inline static void bitmask_set_bit(gmx_bitmask_t* m, int b)
137 (*m)[b / 64] |= (static_cast<uint64_t>(1) << (b % 64));
140 inline static void bitmask_init_bit(gmx_bitmask_t* m, int b)
142 bitmask_clear(m);
143 (*m)[b / 64] = (static_cast<uint64_t>(1) << (b % 64));
146 inline static void bitmask_init_low_bits(gmx_bitmask_t* m, int b)
148 memset(m->data(), 255, b / 64 * 8);
149 (*m)[b / 64] = (static_cast<uint64_t>(1) << (b % 64)) - 1;
150 memset(m->data() + (b / 64 + 1), 0, (BITMASK_ALEN - b / 64 - 1) * 8);
153 inline static bool bitmask_is_set(gmx_bitmask_t m, int b)
155 return (m[b / 64] & (static_cast<uint64_t>(1) << (b % 64))) != 0;
158 inline static bool bitmask_is_disjoint(gmx_bitmask_t a, gmx_bitmask_t b)
160 int i;
161 bool r = true;
162 for (i = 0; i < BITMASK_ALEN; i++)
164 r = r && ((a[i] & b[i]) == 0U);
166 return r;
169 inline static bool bitmask_is_equal(gmx_bitmask_t a, gmx_bitmask_t b)
171 int i;
172 bool r = true;
173 for (i = 0; i < BITMASK_ALEN; i++)
175 r = r && (a[i] == b[i]);
177 return r;
180 inline static bool bitmask_is_zero(gmx_bitmask_t m)
182 int i;
183 bool r = true;
184 for (i = 0; i < BITMASK_ALEN; i++)
186 r = r && (m[i] == 0U);
188 return r;
191 inline static void bitmask_union(gmx_bitmask_t* a, gmx_bitmask_t b)
193 int i;
194 for (i = 0; i < BITMASK_ALEN; i++)
196 (*a)[i] |= b[i];
199 #endif
200 // In bitmask.h because only current use is for bitmask.
202 //! Convert uint32_t to hex string
203 inline static std::string to_hex_string(uint32_t m)
205 return gmx::formatString("%08" PRIx32, m);
207 //! Convert uint64_t to hex string
208 inline static std::string to_hex_string(uint64_t m)
210 return gmx::formatString("%016" PRIx64, m);
212 //! Convert container of intergers to hex string
213 template<typename C>
214 inline static std::string to_hex_string(C m)
216 std::string ret;
217 for (auto it = m.rbegin(); it < m.rend(); ++it)
219 ret += to_hex_string(*it);
221 return ret;
224 #endif