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
44 #include "gromacs/utility/basedefinitions.h"
45 #include "gromacs/utility/gmxassert.h"
46 #include "gromacs/utility/range.h"
51 /*! \brief Division of a range of indices into consecutive blocks
53 * A range of consecutive indices 0 to full.range.end() is divided
54 * into numBlocks() consecutive blocks of consecutive indices.
55 * Block b contains indices i for which block(b).begin() <= i < block(b).end().
57 class RangePartitioning
60 /*! \brief A block defined by a range of atom indices */
61 using Block
= Range
<int>;
63 /*! \brief Returns the number of blocks */
66 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
+ 1]);
75 /*! \brief Returns the full range */
76 Block
fullRange() const
78 return Block(index_
.front(), index_
.back());
81 /*! \brief Returns a range starting at \p blockIndexBegin and ending at \p blockIndexEnd */
82 Block
subRange(int blockIndexBegin
,
83 int blockIndexEnd
) const
85 return Block(index_
[blockIndexBegin
], index_
[blockIndexEnd
]);
88 /*! \brief Returns true when all blocks have size 0 or numBlocks()=0 */
89 bool allBlocksHaveSizeOne() const
91 return (index_
.back() == numBlocks());
94 /*! \brief Appends a block of size \p blockSize at the end of the range
96 * \note blocksize has to be >= 1
98 void appendBlock(int blockSize
)
100 GMX_ASSERT(blockSize
> 0, "block sizes should be >= 1");
101 index_
.push_back(index_
.back() + blockSize
);
104 /*! \brief Removes all blocks */
110 /*! \brief Reduces the number of blocks to \p newNumBlocks
112 * \note \p newNumBlocks should be <= numBlocks().
114 void reduceNumBlocks(int newNumBlocks
)
116 GMX_ASSERT(newNumBlocks
<= numBlocks(), "Can only shrink to fewer blocks");
117 index_
.resize(newNumBlocks
+ 1);
120 /*! \brief Sets the partitioning to \p numBlocks blocks each of size 1 */
121 void setAllBlocksSizeOne(int numBlocks
);
123 /*! \brief Returns the raw block index array, avoid using this */
124 std::vector
<int> &rawIndex()
130 std::vector
<int> index_
= { 0 }; /**< The list of block begin/end indices */
135 /* Deprecated, C-style version of RangePartitioning */
136 typedef struct t_block
138 int blockSize(int blockIndex
) const
140 GMX_ASSERT(blockIndex
< nr
, "blockIndex should be in range");
141 return index
[blockIndex
+ 1] - index
[blockIndex
];
144 int nr
; /* The number of blocks */
145 int *index
; /* Array of indices (dim: nr+1) */
146 int nalloc_index
; /* The allocation size for index */
151 int nr
; /* The number of blocks */
152 int *index
; /* Array of indices in a (dim: nr+1) */
153 int nra
; /* The number of atoms */
154 int *a
; /* Array of atom numbers in each group */
156 /* Block i (0<=i<nr) runs from */
157 /* index[i] to index[i+1]-1. There will */
158 /* allways be an extra entry in index */
159 /* to terminate the table */
160 int nalloc_index
; /* The allocation size for index */
161 int nalloc_a
; /* The allocation size for a */
165 * Fully initialize t_block datastructure.
167 * Initializes a \p block and sets up the first index to zero.
169 * \param[in,out] block datastructure to initialize.
171 void init_block(t_block
*block
);
174 * Fully initialize t_blocka datastructure.
176 * Initializes a \p block and sets up the first index to zero.
177 * The atom number array is initialized to nullptr.
179 * \param[in,out] block datastructure to initialize.
181 void init_blocka(t_blocka
*block
);
184 * In general all t_block datastructures should be avoided
185 * in favour of RangePartitioning. This here is a simple cludge
186 * to use more modern initialization while we move to the use
187 * of RangePartitioning.
191 * Minimal initialization of t_block datastructure.
193 * Performs the equivalent to a snew on a t_block, setting all
194 * values to zero or nullptr. Needed for some cases where the topology
195 * handling expects a block to be valid initialized (e.g. during domain
196 * decomposition) but without the first block set to zero.
198 * \param[in,out] block datastructure to initialize.
200 void init_block_null(t_block
*block
);
203 * Minimal initialization of t_blocka datastructure.
205 * Performs the equivalent to a snew on a t_blocka, setting all
206 * values to zero or nullptr. Needed for some cases where the topology
207 * handling expects a block to be valid initialized (e.g. during domain
208 * decomposition) but without the first block set to zero.
210 * \param[in,out] block datastructure to initialize.
212 void init_blocka_null(t_blocka
*block
);
214 t_blocka
*new_blocka();
215 /* allocate new block */
217 void done_block(t_block
*block
);
218 void done_blocka(t_blocka
*block
);
220 void copy_blocka(const t_blocka
*src
, t_blocka
*dest
);
222 void copy_block(const t_block
*src
, t_block
*dst
);
224 void stupid_fill_block(t_block
*grp
, int natom
, gmx_bool bOneIndexGroup
);
225 /* Fill a block structure with numbers identical to the index
226 * (0, 1, 2, .. natom-1)
227 * If bOneIndexGroup, then all atoms are lumped in one index group,
228 * otherwise there is one atom per index entry
231 void stupid_fill_blocka(t_blocka
*grp
, int natom
);
232 /* Fill a block structure with numbers identical to the index
233 * (0, 1, 2, .. natom-1)
234 * There is one atom per index entry
237 void pr_block(FILE *fp
, int indent
, const char *title
, const t_block
*block
, gmx_bool bShowNumbers
);
238 void pr_blocka(FILE *fp
, int indent
, const char *title
, const t_blocka
*block
, gmx_bool bShowNumbers
);