Use proper doxygen tags in modular simulator
[gromacs.git] / src / gromacs / topology / block.h
blob10fee47a24b40f31275487be8f752c85de506af9
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,2020, 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"
46 #include "gromacs/utility/range.h"
48 namespace gmx
51 template<typename>
52 class ListOfLists;
54 /*! \brief Division of a range of indices into consecutive blocks
56 * A range of consecutive indices 0 to full.range.end() is divided
57 * into numBlocks() consecutive blocks of consecutive indices.
58 * Block b contains indices i for which block(b).begin() <= i < block(b).end().
60 class RangePartitioning
62 public:
63 /*! \brief A block defined by a range of atom indices */
64 using Block = Range<int>;
66 /*! \brief Returns the number of blocks */
67 int numBlocks() const { return static_cast<int>(index_.size()) - 1; }
69 /*! \brief Returns the size of the block with index \p blockIndex */
70 Block block(int blockIndex) const
72 return Block(index_[blockIndex], index_[blockIndex + 1LL]);
75 /*! \brief Returns the full range */
76 Block fullRange() const { return Block(index_.front(), index_.back()); }
78 /*! \brief Returns a range starting at \p blockIndexBegin and ending at \p blockIndexEnd */
79 Block subRange(int blockIndexBegin, int blockIndexEnd) const
81 return Block(index_[blockIndexBegin], index_[blockIndexEnd]);
84 /*! \brief Returns true when all blocks have size 0 or numBlocks()=0 */
85 bool allBlocksHaveSizeOne() const { return (index_.back() == numBlocks()); }
87 /*! \brief Appends a block of size \p blockSize at the end of the range
89 * \note blocksize has to be >= 1
91 void appendBlock(int blockSize)
93 GMX_ASSERT(blockSize > 0, "block sizes should be >= 1");
94 index_.push_back(index_.back() + blockSize);
97 /*! \brief Removes all blocks */
98 void clear() { index_.resize(1); }
100 /*! \brief Reduces the number of blocks to \p newNumBlocks
102 * \note \p newNumBlocks should be <= numBlocks().
104 void reduceNumBlocks(int newNumBlocks)
106 GMX_ASSERT(newNumBlocks <= numBlocks(), "Can only shrink to fewer blocks");
107 index_.resize(newNumBlocks + 1LL);
110 /*! \brief Sets the partitioning to \p numBlocks blocks each of size 1 */
111 void setAllBlocksSizeOne(int numBlocks);
113 /*! \brief Returns the raw block index array, avoid using this */
114 std::vector<int>& rawIndex() { return index_; }
116 private:
117 std::vector<int> index_ = { 0 }; /**< The list of block begin/end indices */
120 } // namespace gmx
122 /* Deprecated, C-style version of RangePartitioning */
123 typedef struct t_block
125 int blockSize(int blockIndex) const
127 GMX_ASSERT(blockIndex < nr, "blockIndex should be in range");
128 return index[blockIndex + 1] - index[blockIndex];
131 int nr; /* The number of blocks */
132 int* index; /* Array of indices (dim: nr+1) */
133 int nalloc_index; /* The allocation size for index */
134 } t_block;
136 struct t_blocka
138 int nr; /* The number of blocks */
139 int* index; /* Array of indices in a (dim: nr+1) */
140 int nra; /* The number of atoms */
141 int* a; /* Array of atom numbers in each group */
142 /* (dim: nra) */
143 /* Block i (0<=i<nr) runs from */
144 /* index[i] to index[i+1]-1. There will */
145 /* allways be an extra entry in index */
146 /* to terminate the table */
147 int nalloc_index; /* The allocation size for index */
148 int nalloc_a; /* The allocation size for a */
151 /*! \brief
152 * Fully initialize t_block datastructure.
154 * Initializes a \p block and sets up the first index to zero.
156 * \param[in,out] block datastructure to initialize.
158 void init_block(t_block* block);
160 /*! \brief
161 * Fully initialize t_blocka datastructure.
163 * Initializes a \p block and sets up the first index to zero.
164 * The atom number array is initialized to nullptr.
166 * \param[in,out] block datastructure to initialize.
168 void init_blocka(t_blocka* block);
170 /* TODO
171 * In general all t_block datastructures should be avoided
172 * in favour of RangePartitioning. This here is a simple cludge
173 * to use more modern initialization while we move to the use
174 * of RangePartitioning.
177 /*! \brief
178 * Minimal initialization of t_block datastructure.
180 * Performs the equivalent to a snew on a t_block, setting all
181 * values to zero or nullptr. Needed for some cases where the topology
182 * handling expects a block to be valid initialized (e.g. during domain
183 * decomposition) but without the first block set to zero.
185 * \param[in,out] block datastructure to initialize.
187 void init_block_null(t_block* block);
189 /*! \brief
190 * Minimal initialization of t_blocka datastructure.
192 * Performs the equivalent to a snew on a t_blocka, setting all
193 * values to zero or nullptr. Needed for some cases where the topology
194 * handling expects a block to be valid initialized (e.g. during domain
195 * decomposition) but without the first block set to zero.
197 * \param[in,out] block datastructure to initialize.
199 void init_blocka_null(t_blocka* block);
201 t_blocka* new_blocka();
202 /* allocate new block */
204 void done_block(t_block* block);
205 void done_blocka(t_blocka* block);
207 void copy_blocka(const t_blocka* src, t_blocka* dest);
209 void copy_block(const t_block* src, t_block* dst);
211 void stupid_fill_block(t_block* grp, int natom, gmx_bool bOneIndexGroup);
212 /* Fill a block structure with numbers identical to the index
213 * (0, 1, 2, .. natom-1)
214 * If bOneIndexGroup, then all atoms are lumped in one index group,
215 * otherwise there is one atom per index entry
218 void stupid_fill_blocka(t_blocka* grp, int natom);
219 /* Fill a block structure with numbers identical to the index
220 * (0, 1, 2, .. natom-1)
221 * There is one atom per index entry
224 void pr_block(FILE* fp, int indent, const char* title, const t_block* block, gmx_bool bShowNumbers);
225 void pr_blocka(FILE* fp, int indent, const char* title, const t_blocka* block, gmx_bool bShowNumbers);
226 void pr_listoflists(FILE* fp, int indent, const char* title, const gmx::ListOfLists<int>* block, gmx_bool bShowNumbers);
228 #endif