Update instructions in containers.rst
[gromacs.git] / src / gromacs / topology / block.cpp
blobcc6ad6b191fb593251125fc68a1569d4bff94ae7
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) 2013,2014,2015,2017,2018 by the GROMACS development team.
7 * Copyright (c) 2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 #include "gmxpre.h"
40 #include "block.h"
42 #include <cstdio>
44 #include <algorithm>
46 #include "gromacs/utility/listoflists.h"
47 #include "gromacs/utility/smalloc.h"
48 #include "gromacs/utility/txtdump.h"
50 void gmx::RangePartitioning::setAllBlocksSizeOne(int numBlocksToSet)
52 if (!allBlocksHaveSizeOne())
54 clear();
56 if (numBlocksToSet < numBlocks())
58 index_.resize(numBlocksToSet + 1);
60 else if (numBlocksToSet > numBlocks())
62 for (int b = numBlocks(); b < numBlocksToSet; b++)
64 appendBlock(1);
69 void init_block(t_block* block)
71 block->nr = 0;
72 block->nalloc_index = 1;
73 snew(block->index, block->nalloc_index);
74 block->index[0] = 0;
77 void init_block_null(t_block* block)
79 block->nr = 0;
80 block->nalloc_index = 0;
81 block->index = nullptr;
84 void init_blocka(t_blocka* block)
86 block->nr = 0;
87 block->nra = 0;
88 block->nalloc_index = 1;
89 snew(block->index, block->nalloc_index);
90 block->index[0] = 0;
91 block->nalloc_a = 0;
92 block->a = nullptr;
95 void init_blocka_null(t_blocka* block)
97 block->nr = 0;
98 block->nra = 0;
99 block->nalloc_index = 0;
100 block->index = nullptr;
101 block->nalloc_a = 0;
102 block->a = nullptr;
105 t_blocka* new_blocka()
107 t_blocka* block;
109 snew(block, 1);
110 snew(block->index, 1);
112 return block;
115 void done_block(t_block* block)
117 block->nr = 0;
118 sfree(block->index);
119 block->index = nullptr;
120 block->nalloc_index = 0;
123 void done_blocka(t_blocka* block)
125 block->nr = 0;
126 block->nra = 0;
127 sfree(block->index);
128 sfree(block->a);
129 block->index = nullptr;
130 block->a = nullptr;
131 block->nalloc_index = 0;
132 block->nalloc_a = 0;
135 void stupid_fill_block(t_block* grp, int natom, gmx_bool bOneIndexGroup)
137 if (bOneIndexGroup)
139 grp->nalloc_index = 2;
140 srenew(grp->index, grp->nalloc_index);
141 grp->index[0] = 0;
142 grp->index[1] = natom;
143 grp->nr = 1;
145 else
147 grp->nalloc_index = natom + 1;
148 srenew(grp->index, grp->nalloc_index);
149 for (int i = 0; i <= natom; ++i)
151 grp->index[i] = i;
153 grp->nr = natom;
157 void stupid_fill_blocka(t_blocka* grp, int natom)
159 grp->nalloc_a = natom;
160 snew(grp->a, grp->nalloc_a);
161 for (int i = 0; i < natom; ++i)
163 grp->a[i] = i;
165 grp->nra = natom;
167 grp->nalloc_index = natom + 1;
168 snew(grp->index, grp->nalloc_index);
169 for (int i = 0; i <= natom; ++i)
171 grp->index[i] = i;
173 grp->nr = natom;
176 void copy_blocka(const t_blocka* src, t_blocka* dest)
178 dest->nr = src->nr;
179 /* Workaround for inconsistent handling of nalloc_index in
180 * other parts of the code. Often nalloc_index and nalloc_a
181 * are not set.
183 dest->nalloc_index = std::max(src->nalloc_index, dest->nr + 1);
184 snew(dest->index, dest->nalloc_index);
185 for (int i = 0; i < dest->nr + 1; ++i)
187 dest->index[i] = src->index[i];
189 dest->nra = src->nra;
190 /* See above. */
191 dest->nalloc_a = std::max(src->nalloc_a, dest->nra);
192 snew(dest->a, dest->nalloc_a);
193 for (int i = 0; i < dest->nra; ++i)
195 dest->a[i] = src->a[i];
199 static int pr_block_title(FILE* fp, int indent, const char* title, const t_block* block)
201 if (available(fp, block, indent, title))
203 indent = pr_title(fp, indent, title);
204 pr_indent(fp, indent);
205 fprintf(fp, "nr=%d\n", block->nr);
207 return indent;
210 static int pr_blocka_title(FILE* fp, int indent, const char* title, const t_blocka* block)
212 if (available(fp, block, indent, title))
214 indent = pr_title(fp, indent, title);
215 pr_indent(fp, indent);
216 fprintf(fp, "nr=%d\n", block->nr);
217 pr_indent(fp, indent);
218 fprintf(fp, "nra=%d\n", block->nra);
220 return indent;
223 static int pr_listoflists_title(FILE* fp, int indent, const char* title, const gmx::ListOfLists<int>* lists)
225 if (available(fp, lists, indent, title))
227 indent = pr_title(fp, indent, title);
228 pr_indent(fp, indent);
229 fprintf(fp, "numLists=%zu\n", lists->size());
230 pr_indent(fp, indent);
231 fprintf(fp, "numElements=%d\n", lists->numElements());
233 return indent;
236 static void low_pr_blocka(FILE* fp, int indent, const char* title, const t_blocka* block, gmx_bool bShowNumbers)
238 int i;
240 if (available(fp, block, indent, title))
242 indent = pr_blocka_title(fp, indent, title, block);
243 for (i = 0; i <= block->nr; i++)
245 pr_indent(fp, indent + INDENT);
246 fprintf(fp, "%s->index[%d]=%d\n", title, bShowNumbers ? i : -1, block->index[i]);
248 for (i = 0; i < block->nra; i++)
250 pr_indent(fp, indent + INDENT);
251 fprintf(fp, "%s->a[%d]=%d\n", title, bShowNumbers ? i : -1, block->a[i]);
256 void pr_block(FILE* fp, int indent, const char* title, const t_block* block, gmx_bool bShowNumbers)
258 int i, start;
260 if (available(fp, block, indent, title))
262 indent = pr_block_title(fp, indent, title, block);
263 start = 0;
264 if (block->index[start] != 0)
266 fprintf(fp, "block->index[%d] should be 0\n", start);
268 else
270 for (i = 0; i < block->nr; i++)
272 int end = block->index[i + 1];
273 pr_indent(fp, indent);
274 if (end <= start)
276 fprintf(fp, "%s[%d]={}\n", title, i);
278 else
280 fprintf(fp, "%s[%d]={%d..%d}\n", title, bShowNumbers ? i : -1,
281 bShowNumbers ? start : -1, bShowNumbers ? end - 1 : -1);
283 start = end;
289 void pr_blocka(FILE* fp, int indent, const char* title, const t_blocka* block, gmx_bool bShowNumbers)
291 int i, j, ok, size, start, end;
293 if (available(fp, block, indent, title))
295 indent = pr_blocka_title(fp, indent, title, block);
296 start = 0;
297 end = start;
298 if ((ok = static_cast<int>(block->index[start] == 0)) == 0)
300 fprintf(fp, "block->index[%d] should be 0\n", start);
302 else
304 for (i = 0; i < block->nr; i++)
306 end = block->index[i + 1];
307 size = pr_indent(fp, indent);
308 if (end <= start)
310 size += fprintf(fp, "%s[%d]={", title, i);
312 else
314 size += fprintf(fp, "%s[%d][%d..%d]={", title, bShowNumbers ? i : -1,
315 bShowNumbers ? start : -1, bShowNumbers ? end - 1 : -1);
317 for (j = start; j < end; j++)
319 if (j > start)
321 size += fprintf(fp, ", ");
323 if ((size) > (USE_WIDTH))
325 fprintf(fp, "\n");
326 size = pr_indent(fp, indent + INDENT);
328 size += fprintf(fp, "%d", block->a[j]);
330 fprintf(fp, "}\n");
331 start = end;
334 if ((end != block->nra) || (!ok))
336 pr_indent(fp, indent);
337 fprintf(fp, "tables inconsistent, dumping complete tables:\n");
338 low_pr_blocka(fp, indent, title, block, bShowNumbers);
343 void pr_listoflists(FILE* fp, int indent, const char* title, const gmx::ListOfLists<int>* lists, gmx_bool bShowNumbers)
345 if (available(fp, lists, indent, title))
347 indent = pr_listoflists_title(fp, indent, title, lists);
348 for (gmx::index i = 0; i < lists->ssize(); i++)
350 int size = pr_indent(fp, indent);
351 gmx::ArrayRef<const int> list = (*lists)[i];
352 if (list.empty())
354 size += fprintf(fp, "%s[%d]={", title, int(i));
356 else
358 size += fprintf(fp, "%s[%d][num=%zu]={", title, bShowNumbers ? int(i) : -1, list.size());
360 bool isFirst = true;
361 for (const int j : list)
363 if (!isFirst)
365 size += fprintf(fp, ", ");
367 if ((size) > (USE_WIDTH))
369 fprintf(fp, "\n");
370 size = pr_indent(fp, indent + INDENT);
372 size += fprintf(fp, "%d", j);
373 isFirst = false;
375 fprintf(fp, "}\n");
380 void copy_block(const t_block* src, t_block* dst)
382 dst->nr = src->nr;
383 /* Workaround for inconsistent handling of nalloc_index in
384 * other parts of the code. Often nalloc_index and nalloc_a
385 * are not set.
387 dst->nalloc_index = std::max(src->nalloc_index, dst->nr + 1);
388 snew(dst->index, dst->nalloc_index);
389 for (int i = 0; i < dst->nr + 1; ++i)
391 dst->index[i] = src->index[i];