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) 2013,2014,2015,2017,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.
45 #include "gromacs/utility/smalloc.h"
46 #include "gromacs/utility/txtdump.h"
48 void gmx::RangePartitioning::setAllBlocksSizeOne(int numBlocksToSet
)
50 if (!allBlocksHaveSizeOne())
54 if (numBlocksToSet
< numBlocks())
56 index_
.resize(numBlocksToSet
+ 1);
58 else if (numBlocksToSet
> numBlocks())
60 for (int b
= numBlocks(); b
< numBlocksToSet
; b
++)
67 void init_block(t_block
*block
)
70 block
->nalloc_index
= 1;
71 snew(block
->index
, block
->nalloc_index
);
75 void init_blocka(t_blocka
*block
)
79 block
->nalloc_index
= 1;
80 snew(block
->index
, block
->nalloc_index
);
86 t_blocka
*new_blocka()
91 snew(block
->index
, 1);
96 void done_block(t_block
*block
)
100 block
->index
= nullptr;
101 block
->nalloc_index
= 0;
104 void done_blocka(t_blocka
*block
)
110 block
->index
= nullptr;
112 block
->nalloc_index
= 0;
116 void stupid_fill_block(t_block
*grp
, int natom
, gmx_bool bOneIndexGroup
)
120 grp
->nalloc_index
= 2;
121 srenew(grp
->index
, grp
->nalloc_index
);
123 grp
->index
[1] = natom
;
128 grp
->nalloc_index
= natom
+1;
129 srenew(grp
->index
, grp
->nalloc_index
);
130 for (int i
= 0; i
<= natom
; ++i
)
138 void stupid_fill_blocka(t_blocka
*grp
, int natom
)
140 grp
->nalloc_a
= natom
;
141 snew(grp
->a
, grp
->nalloc_a
);
142 for (int i
= 0; i
< natom
; ++i
)
148 grp
->nalloc_index
= natom
+ 1;
149 snew(grp
->index
, grp
->nalloc_index
);
150 for (int i
= 0; i
<= natom
; ++i
)
157 void copy_blocka(const t_blocka
*src
, t_blocka
*dest
)
160 /* Workaround for inconsistent handling of nalloc_index in
161 * other parts of the code. Often nalloc_index and nalloc_a
164 dest
->nalloc_index
= std::max(src
->nalloc_index
, dest
->nr
+ 1);
165 snew(dest
->index
, dest
->nalloc_index
);
166 for (int i
= 0; i
< dest
->nr
+1; ++i
)
168 dest
->index
[i
] = src
->index
[i
];
170 dest
->nra
= src
->nra
;
172 dest
->nalloc_a
= std::max(src
->nalloc_a
, dest
->nra
);
173 snew(dest
->a
, dest
->nalloc_a
);
174 for (int i
= 0; i
< dest
->nra
; ++i
)
176 dest
->a
[i
] = src
->a
[i
];
180 static int pr_block_title(FILE *fp
, int indent
, const char *title
, const t_block
*block
)
182 if (available(fp
, block
, indent
, title
))
184 indent
= pr_title(fp
, indent
, title
);
185 pr_indent(fp
, indent
);
186 fprintf(fp
, "nr=%d\n", block
->nr
);
191 static int pr_blocka_title(FILE *fp
, int indent
, const char *title
, const t_blocka
*block
)
193 if (available(fp
, block
, indent
, title
))
195 indent
= pr_title(fp
, indent
, title
);
196 pr_indent(fp
, indent
);
197 fprintf(fp
, "nr=%d\n", block
->nr
);
198 pr_indent(fp
, indent
);
199 fprintf(fp
, "nra=%d\n", block
->nra
);
204 static void low_pr_blocka(FILE *fp
, int indent
, const char *title
, const t_blocka
*block
, gmx_bool bShowNumbers
)
208 if (available(fp
, block
, indent
, title
))
210 indent
= pr_blocka_title(fp
, indent
, title
, block
);
211 for (i
= 0; i
<= block
->nr
; i
++)
213 pr_indent(fp
, indent
+INDENT
);
214 fprintf(fp
, "%s->index[%d]=%d\n",
215 title
, bShowNumbers
? i
: -1, block
->index
[i
]);
217 for (i
= 0; i
< block
->nra
; i
++)
219 pr_indent(fp
, indent
+INDENT
);
220 fprintf(fp
, "%s->a[%d]=%d\n",
221 title
, bShowNumbers
? i
: -1, block
->a
[i
]);
226 void pr_block(FILE *fp
, int indent
, const char *title
, const t_block
*block
, gmx_bool bShowNumbers
)
230 if (available(fp
, block
, indent
, title
))
232 indent
= pr_block_title(fp
, indent
, title
, block
);
234 if (block
->index
[start
] != 0)
236 fprintf(fp
, "block->index[%d] should be 0\n", start
);
240 for (i
= 0; i
< block
->nr
; i
++)
242 int end
= block
->index
[i
+1];
243 pr_indent(fp
, indent
);
246 fprintf(fp
, "%s[%d]={}\n", title
, i
);
250 fprintf(fp
, "%s[%d]={%d..%d}\n",
251 title
, bShowNumbers
? i
: -1,
252 bShowNumbers
? start
: -1, bShowNumbers
? end
-1 : -1);
260 void pr_blocka(FILE *fp
, int indent
, const char *title
, const t_blocka
*block
, gmx_bool bShowNumbers
)
262 int i
, j
, ok
, size
, start
, end
;
264 if (available(fp
, block
, indent
, title
))
266 indent
= pr_blocka_title(fp
, indent
, title
, block
);
269 if ((ok
= static_cast<int>(block
->index
[start
] == 0)) == 0)
271 fprintf(fp
, "block->index[%d] should be 0\n", start
);
275 for (i
= 0; i
< block
->nr
; i
++)
277 end
= block
->index
[i
+1];
278 size
= pr_indent(fp
, indent
);
281 size
+= fprintf(fp
, "%s[%d]={", title
, i
);
285 size
+= fprintf(fp
, "%s[%d][%d..%d]={",
286 title
, bShowNumbers
? i
: -1,
287 bShowNumbers
? start
: -1, bShowNumbers
? end
-1 : -1);
289 for (j
= start
; j
< end
; j
++)
293 size
+= fprintf(fp
, ", ");
295 if ((size
) > (USE_WIDTH
))
298 size
= pr_indent(fp
, indent
+INDENT
);
300 size
+= fprintf(fp
, "%d", block
->a
[j
]);
306 if ((end
!= block
->nra
) || (!ok
))
308 pr_indent(fp
, indent
);
309 fprintf(fp
, "tables inconsistent, dumping complete tables:\n");
310 low_pr_blocka(fp
, indent
, title
, block
, bShowNumbers
);