Update clang-tidy to clang version 8
[gromacs.git] / src / gromacs / topology / block.h
blobe39e0f72f4d0aea6b19241b372cfa9b00880890e
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2010,2014,2015,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 #ifndef GMX_TOPOLOGY_BLOCK_H
38 #define GMX_TOPOLOGY_BLOCK_H
40 #include <stdio.h>
42 #include <vector>
44 #include "gromacs/utility/basedefinitions.h"
45 #include "gromacs/utility/gmxassert.h"
47 namespace gmx
50 /*! \brief Division of a range of indices into consecutive blocks
52 * A range of consecutive indices 0 to full.range.end() is divided
53 * into numBlocks() consecutive blocks of consecutive indices.
54 * Block b contains indices i for which block(b).begin() <= i < block(b).end().
56 class RangePartitioning
58 public:
59 /*! \brief Struct for returning the range of a block.
61 * Can be used in a range loop.
63 struct Block
65 public:
66 /*! \brief An iterator that loops over integers */
67 struct iterator
69 //! Constructor
70 iterator(int value) : value_(value) {}
71 //! Value
72 operator int () const { return value_; }
73 //! Reference
74 operator int &() { return value_; }
75 //! Pointer
76 int operator* () const { return value_; }
77 //! Inequality comparison
78 bool operator!= (const iterator other) { return value_ != other; }
79 //! Increment operator
80 iterator &operator++() { ++value_; return *this; }
81 //! Increment operator
82 iterator operator++(int) { iterator tmp(*this); ++value_; return tmp; }
83 //! The actual value
84 int value_;
87 /*! \brief Constructor, constructs a range starting at 0 with 0 blocks */
88 Block(int begin,
89 int end) :
90 begin_(begin),
91 end_(end)
95 /*! \brief Begin iterator/value */
96 iterator begin() const { return begin_; }
97 /*! \brief End iterator/value */
98 iterator end() const { return end_; }
100 /*! \brief The number of items in the block */
101 int size() const
103 return end_ - begin_;
106 /*! \brief Returns whether \p index is within range of the block */
107 bool inRange(int index) const
109 return (begin_ <= index && index < end_);
112 private:
113 const int begin_; /**< The start index of the block */
114 const int end_; /**< The end index of the block */
117 /*! \brief Returns the number of blocks */
118 int numBlocks() const
120 return static_cast<int>(index_.size()) - 1;
123 /*! \brief Returns the size of the block with index \p blockIndex */
124 Block block(int blockIndex) const
126 return Block(index_[blockIndex], index_[blockIndex + 1]);
129 /*! \brief Returns the full range */
130 Block fullRange() const
132 return Block(index_.front(), index_.back());
135 /*! \brief Returns a range starting at \p blockIndexBegin and ending at \p blockIndexEnd */
136 Block subRange(int blockIndexBegin,
137 int blockIndexEnd) const
139 return Block(index_[blockIndexBegin], index_[blockIndexEnd]);
142 /*! \brief Returns true when all blocks have size 0 or numBlocks()=0 */
143 bool allBlocksHaveSizeOne() const
145 return (index_.back() == numBlocks());
148 /*! \brief Appends a block of size \p blockSize at the end of the range
150 * \note blocksize has to be >= 1
152 void appendBlock(int blockSize)
154 GMX_ASSERT(blockSize > 0, "block sizes should be >= 1");
155 index_.push_back(index_.back() + blockSize);
158 /*! \brief Removes all blocks */
159 void clear()
161 index_.resize(1);
164 /*! \brief Reduces the number of blocks to \p newNumBlocks
166 * \note \p newNumBlocks should be <= numBlocks().
168 void reduceNumBlocks(int newNumBlocks)
170 GMX_ASSERT(newNumBlocks <= numBlocks(), "Can only shrink to fewer blocks");
171 index_.resize(newNumBlocks + 1);
174 /*! \brief Sets the partitioning to \p numBlocks blocks each of size 1 */
175 void setAllBlocksSizeOne(int numBlocks);
177 /*! \brief Returns the raw block index array, avoid using this */
178 std::vector<int> &rawIndex()
180 return index_;
183 private:
184 std::vector<int> index_ = { 0 }; /**< The list of block begin/end indices */
187 } // namespace gmx
189 /* Deprecated, C-style version of RangePartitioning */
190 typedef struct t_block
192 int blockSize(int blockIndex) const
194 GMX_ASSERT(blockIndex < nr, "blockIndex should be in range");
195 return index[blockIndex + 1] - index[blockIndex];
198 int nr; /* The number of blocks */
199 int *index; /* Array of indices (dim: nr+1) */
200 int nalloc_index; /* The allocation size for index */
201 } t_block;
203 struct t_blocka
205 int nr; /* The number of blocks */
206 int *index; /* Array of indices in a (dim: nr+1) */
207 int nra; /* The number of atoms */
208 int *a; /* Array of atom numbers in each group */
209 /* (dim: nra) */
210 /* Block i (0<=i<nr) runs from */
211 /* index[i] to index[i+1]-1. There will */
212 /* allways be an extra entry in index */
213 /* to terminate the table */
214 int nalloc_index; /* The allocation size for index */
215 int nalloc_a; /* The allocation size for a */
218 /*! \brief
219 * Fully initialize t_block datastructure.
221 * Initializes a \p block and sets up the first index to zero.
223 * \param[in,out] block datastructure to initialize.
225 void init_block(t_block *block);
227 /*! \brief
228 * Fully initialize t_blocka datastructure.
230 * Initializes a \p block and sets up the first index to zero.
231 * The atom number array is initialized to nullptr.
233 * \param[in,out] block datastructure to initialize.
235 void init_blocka(t_blocka *block);
237 /* TODO
238 * In general all t_block datastructures should be avoided
239 * in favour of RangePartitioning. This here is a simple cludge
240 * to use more modern initialization while we move to the use
241 * of RangePartitioning.
244 /*! \brief
245 * Minimal initialization of t_block datastructure.
247 * Performs the equivalent to a snew on a t_block, setting all
248 * values to zero or nullptr. Needed for some cases where the topology
249 * handling expects a block to be valid initialized (e.g. during domain
250 * decomposition) but without the first block set to zero.
252 * \param[in,out] block datastructure to initialize.
254 void init_block_null(t_block *block);
256 /*! \brief
257 * Minimal initialization of t_blocka datastructure.
259 * Performs the equivalent to a snew on a t_blocka, setting all
260 * values to zero or nullptr. Needed for some cases where the topology
261 * handling expects a block to be valid initialized (e.g. during domain
262 * decomposition) but without the first block set to zero.
264 * \param[in,out] block datastructure to initialize.
266 void init_blocka_null(t_blocka *block);
268 t_blocka *new_blocka();
269 /* allocate new block */
271 void done_block(t_block *block);
272 void done_blocka(t_blocka *block);
274 void copy_blocka(const t_blocka *src, t_blocka *dest);
276 void copy_block(const t_block *src, t_block *dst);
278 void stupid_fill_block(t_block *grp, int natom, gmx_bool bOneIndexGroup);
279 /* Fill a block structure with numbers identical to the index
280 * (0, 1, 2, .. natom-1)
281 * If bOneIndexGroup, then all atoms are lumped in one index group,
282 * otherwise there is one atom per index entry
285 void stupid_fill_blocka(t_blocka *grp, int natom);
286 /* Fill a block structure with numbers identical to the index
287 * (0, 1, 2, .. natom-1)
288 * There is one atom per index entry
291 void pr_block(FILE *fp, int indent, const char *title, const t_block *block, gmx_bool bShowNumbers);
292 void pr_blocka(FILE *fp, int indent, const char *title, const t_blocka *block, gmx_bool bShowNumbers);
294 #endif