Remove unused thole polarization rfac parameter
[alexxy/gromacs.git] / src / gromacs / topology / block.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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,2021, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 #ifndef GMX_TOPOLOGY_BLOCK_H
38 #define GMX_TOPOLOGY_BLOCK_H
39
40 #include <stdio.h>
41
42 #include <vector>
43
44 #include "gromacs/utility/gmxassert.h"
45 #include "gromacs/utility/range.h"
46
47 namespace gmx
48 {
49
50 template<typename>
51 class ListOfLists;
52
53 /*! \brief Division of a range of indices into consecutive blocks
54  *
55  * A range of consecutive indices 0 to full.range.end() is divided
56  * into numBlocks() consecutive blocks of consecutive indices.
57  * Block b contains indices i for which block(b).begin() <= i < block(b).end().
58  */
59 class RangePartitioning
60 {
61 public:
62     /*! \brief A block defined by a range of atom indices */
63     using Block = Range<int>;
64
65     /*! \brief Returns the number of blocks */
66     int numBlocks() const { return static_cast<int>(index_.size()) - 1; }
67
68     /*! \brief Returns the size of the block with index \p blockIndex */
69     Block block(int blockIndex) const
70     {
71         return Block(index_[blockIndex], index_[blockIndex + 1LL]);
72     }
73
74     /*! \brief Returns the full range */
75     Block fullRange() const { return Block(index_.front(), index_.back()); }
76
77     /*! \brief Returns a range starting at \p blockIndexBegin and ending at \p blockIndexEnd */
78     Block subRange(int blockIndexBegin, int blockIndexEnd) const
79     {
80         return Block(index_[blockIndexBegin], index_[blockIndexEnd]);
81     }
82
83     /*! \brief Returns true when all blocks have size 0 or numBlocks()=0 */
84     bool allBlocksHaveSizeOne() const { return (index_.back() == numBlocks()); }
85
86     /*! \brief Appends a block of size \p blockSize at the end of the range
87      *
88      * \note blocksize has to be >= 1
89      */
90     void appendBlock(int blockSize)
91     {
92         GMX_ASSERT(blockSize > 0, "block sizes should be >= 1");
93         index_.push_back(index_.back() + blockSize);
94     }
95
96     /*! \brief Removes all blocks */
97     void clear() { index_.resize(1); }
98
99     /*! \brief Reduces the number of blocks to \p newNumBlocks
100      *
101      * \note \p newNumBlocks should be <= numBlocks().
102      */
103     void reduceNumBlocks(int newNumBlocks)
104     {
105         GMX_ASSERT(newNumBlocks <= numBlocks(), "Can only shrink to fewer blocks");
106         index_.resize(newNumBlocks + 1LL);
107     }
108
109     /*! \brief Sets the partitioning to \p numBlocks blocks each of size 1 */
110     void setAllBlocksSizeOne(int numBlocks);
111
112     /*! \brief Returns the raw block index array, avoid using this */
113     std::vector<int>& rawIndex() { return index_; }
114
115 private:
116     std::vector<int> index_ = { 0 }; /**< The list of block begin/end indices */
117 };
118
119 } // namespace gmx
120
121 /* Deprecated, C-style version of RangePartitioning */
122 typedef struct t_block
123 {
124     int blockSize(int blockIndex) const
125     {
126         GMX_ASSERT(blockIndex < nr, "blockIndex should be in range");
127         return index[blockIndex + 1] - index[blockIndex];
128     }
129
130     int  nr;           /* The number of blocks          */
131     int* index;        /* Array of indices (dim: nr+1)  */
132     int  nalloc_index; /* The allocation size for index */
133 } t_block;
134
135 struct t_blocka
136 {
137     int  nr;    /* The number of blocks              */
138     int* index; /* Array of indices in a (dim: nr+1) */
139     int  nra;   /* The number of atoms               */
140     int* a;     /* Array of atom numbers in each group  */
141     /* (dim: nra)                           */
142     /* Block i (0<=i<nr) runs from          */
143     /* index[i] to index[i+1]-1. There will */
144     /* allways be an extra entry in index   */
145     /* to terminate the table               */
146     int nalloc_index; /* The allocation size for index        */
147     int nalloc_a;     /* The allocation size for a            */
148 };
149
150 /*! \brief
151  * Fully initialize t_block datastructure.
152  *
153  * Initializes a \p block and sets up the first index to zero.
154  *
155  * \param[in,out] block datastructure to initialize.
156  */
157 void init_block(t_block* block);
158
159 /*! \brief
160  * Fully initialize t_blocka datastructure.
161  *
162  * Initializes a \p block and sets up the first index to zero.
163  * The atom number array is initialized to nullptr.
164  *
165  * \param[in,out] block datastructure to initialize.
166  */
167 void init_blocka(t_blocka* block);
168
169 /* TODO
170  * In general all t_block datastructures should be avoided
171  * in favour of RangePartitioning. This here is a simple cludge
172  * to use more modern initialization while we move to the use
173  * of RangePartitioning.
174  */
175
176 /*! \brief
177  * Minimal initialization of t_block datastructure.
178  *
179  * Performs the equivalent to a snew on a t_block, setting all
180  * values to zero or nullptr. Needed for some cases where the topology
181  * handling expects a block to be valid initialized (e.g. during domain
182  * decomposition) but without the first block set to zero.
183  *
184  * \param[in,out] block datastructure to initialize.
185  */
186 void init_block_null(t_block* block);
187
188 /*! \brief
189  * Minimal initialization of t_blocka datastructure.
190  *
191  * Performs the equivalent to a snew on a t_blocka, setting all
192  * values to zero or nullptr. Needed for some cases where the topology
193  * handling expects a block to be valid initialized (e.g. during domain
194  * decomposition) but without the first block set to zero.
195  *
196  * \param[in,out] block datastructure to initialize.
197  */
198 void init_blocka_null(t_blocka* block);
199
200 t_blocka* new_blocka();
201 /* allocate new block */
202
203 void done_block(t_block* block);
204 void done_blocka(t_blocka* block);
205
206 void copy_blocka(const t_blocka* src, t_blocka* dest);
207
208 void copy_block(const t_block* src, t_block* dst);
209
210 void stupid_fill_block(t_block* grp, int natom, gmx_bool bOneIndexGroup);
211 /* Fill a block structure with numbers identical to the index
212  * (0, 1, 2, .. natom-1)
213  * If bOneIndexGroup, then all atoms are  lumped in one index group,
214  * otherwise there is one atom per index entry
215  */
216
217 void stupid_fill_blocka(t_blocka* grp, int natom);
218 /* Fill a block structure with numbers identical to the index
219  * (0, 1, 2, .. natom-1)
220  * There is one atom per index entry
221  */
222
223 void pr_block(FILE* fp, int indent, const char* title, const t_block* block, gmx_bool bShowNumbers);
224 void pr_blocka(FILE* fp, int indent, const char* title, const t_blocka* block, gmx_bool bShowNumbers);
225 void pr_listoflists(FILE* fp, int indent, const char* title, const gmx::ListOfLists<int>* block, gmx_bool bShowNumbers);
226
227 #endif