Remove mols from gmx_mtop_t
[gromacs.git] / src / gromacs / topology / block.h
blob0131f185ba416aa8fe20131426a19849fbb9d830
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, 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 #ifdef __cplusplus
43 #include <vector>
44 #endif
46 #include "gromacs/utility/basedefinitions.h"
47 #ifdef __cplusplus
48 #include "gromacs/utility/gmxassert.h"
49 #endif
51 #ifdef __cplusplus
52 extern "C" {
53 #endif
55 #ifdef __cplusplus
56 namespace gmx
59 /*! \brief Division of a range of indices into consecutive blocks
61 * A range of consecutive indices 0 to index[numBlocks()] is divided
62 * into numBlocks() consecutive blocks of consecutive indices.
63 * Block b contains indices i for which index[b] <= i < index[b+1].
65 struct BlockRanges
67 /*! \brief Returns the number of blocks
69 * This should only be called on a valid struct. Validy is asserted
70 * (only) in debug mode.
72 int numBlocks() const
74 GMX_ASSERT(index.size() > 0, "numBlocks() should only be called on a valid BlockRanges struct");
75 return index.size() - 1;
78 std::vector<int> index; /**< The list of block begin/end indices */
81 } // nsamespace gmx
82 #endif // __cplusplus
84 /* Deprecated, C-style version of BlockRanges */
85 typedef struct t_block
87 int nr; /* The number of blocks */
88 int *index; /* Array of indices (dim: nr+1) */
89 int nalloc_index; /* The allocation size for index */
90 } t_block;
92 typedef struct t_blocka
94 int nr; /* The number of blocks */
95 int *index; /* Array of indices in a (dim: nr+1) */
96 int nra; /* The number of atoms */
97 int *a; /* Array of atom numbers in each group */
98 /* (dim: nra) */
99 /* Block i (0<=i<nr) runs from */
100 /* index[i] to index[i+1]-1. There will */
101 /* allways be an extra entry in index */
102 /* to terminate the table */
103 int nalloc_index; /* The allocation size for index */
104 int nalloc_a; /* The allocation size for a */
105 } t_blocka;
107 void init_block(t_block *block);
108 void init_blocka(t_blocka *block);
109 t_blocka *new_blocka(void);
110 /* allocate new block */
112 void done_block(t_block *block);
113 void done_blocka(t_blocka *block);
115 void copy_blocka(const t_blocka *src, t_blocka *dest);
117 void stupid_fill_block(t_block *grp, int natom, gmx_bool bOneIndexGroup);
118 /* Fill a block structure with numbers identical to the index
119 * (0, 1, 2, .. natom-1)
120 * If bOneIndexGroup, then all atoms are lumped in one index group,
121 * otherwise there is one atom per index entry
124 void stupid_fill_blocka(t_blocka *grp, int natom);
125 /* Fill a block structure with numbers identical to the index
126 * (0, 1, 2, .. natom-1)
127 * There is one atom per index entry
130 void pr_block(FILE *fp, int indent, const char *title, const t_block *block, gmx_bool bShowNumbers);
131 void pr_blocka(FILE *fp, int indent, const char *title, const t_blocka *block, gmx_bool bShowNumbers);
133 #ifdef __cplusplus
135 #endif
137 #endif