Cleaned up memory usage in gmx traj and trjconv
[gromacs.git] / src / gromacs / topology / block.cpp
blob73c7e47cbbaabb17dec2e52f7ad45eb71081474c
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, 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 #include "gmxpre.h"
39 #include "block.h"
41 #include <cstdio>
43 #include <algorithm>
45 #include "gromacs/utility/smalloc.h"
46 #include "gromacs/utility/txtdump.h"
48 void init_block(t_block *block)
50 block->nr = 0;
51 block->nalloc_index = 1;
52 snew(block->index, block->nalloc_index);
53 block->index[0] = 0;
56 void init_blocka(t_blocka *block)
58 block->nr = 0;
59 block->nra = 0;
60 block->nalloc_index = 1;
61 snew(block->index, block->nalloc_index);
62 block->index[0] = 0;
63 block->nalloc_a = 0;
64 block->a = nullptr;
67 t_blocka *new_blocka(void)
69 t_blocka *block;
71 snew(block, 1);
72 snew(block->index, 1);
74 return block;
77 void done_block(t_block *block)
79 block->nr = 0;
80 sfree(block->index);
81 block->nalloc_index = 0;
84 void done_blocka(t_blocka *block)
86 block->nr = 0;
87 block->nra = 0;
88 sfree(block->index);
89 sfree(block->a);
90 block->index = nullptr;
91 block->a = nullptr;
92 block->nalloc_index = 0;
93 block->nalloc_a = 0;
96 void stupid_fill_block(t_block *grp, int natom, gmx_bool bOneIndexGroup)
98 if (bOneIndexGroup)
100 grp->nalloc_index = 2;
101 snew(grp->index, grp->nalloc_index);
102 grp->index[0] = 0;
103 grp->index[1] = natom;
104 grp->nr = 1;
106 else
108 grp->nalloc_index = natom+1;
109 snew(grp->index, grp->nalloc_index);
110 snew(grp->index, natom+1);
111 for (int i = 0; i <= natom; ++i)
113 grp->index[i] = i;
115 grp->nr = natom;
119 void stupid_fill_blocka(t_blocka *grp, int natom)
121 grp->nalloc_a = natom;
122 snew(grp->a, grp->nalloc_a);
123 for (int i = 0; i < natom; ++i)
125 grp->a[i] = i;
127 grp->nra = natom;
129 grp->nalloc_index = natom + 1;
130 snew(grp->index, grp->nalloc_index);
131 for (int i = 0; i <= natom; ++i)
133 grp->index[i] = i;
135 grp->nr = natom;
138 void copy_blocka(const t_blocka *src, t_blocka *dest)
140 dest->nr = src->nr;
141 /* Workaround for inconsistent handling of nalloc_index in
142 * other parts of the code. Often nalloc_index and nalloc_a
143 * are not set.
145 dest->nalloc_index = std::max(src->nalloc_index, dest->nr + 1);
146 snew(dest->index, dest->nalloc_index);
147 for (int i = 0; i < dest->nr+1; ++i)
149 dest->index[i] = src->index[i];
151 dest->nra = src->nra;
152 /* See above. */
153 dest->nalloc_a = std::max(src->nalloc_a, dest->nra);
154 snew(dest->a, dest->nalloc_a);
155 for (int i = 0; i < dest->nra; ++i)
157 dest->a[i] = src->a[i];
161 static int pr_block_title(FILE *fp, int indent, const char *title, const t_block *block)
163 if (available(fp, block, indent, title))
165 indent = pr_title(fp, indent, title);
166 pr_indent(fp, indent);
167 fprintf(fp, "nr=%d\n", block->nr);
169 return indent;
172 static int pr_blocka_title(FILE *fp, int indent, const char *title, const t_blocka *block)
174 if (available(fp, block, indent, title))
176 indent = pr_title(fp, indent, title);
177 pr_indent(fp, indent);
178 fprintf(fp, "nr=%d\n", block->nr);
179 pr_indent(fp, indent);
180 fprintf(fp, "nra=%d\n", block->nra);
182 return indent;
185 static void low_pr_blocka(FILE *fp, int indent, const char *title, const t_blocka *block, gmx_bool bShowNumbers)
187 int i;
189 if (available(fp, block, indent, title))
191 indent = pr_blocka_title(fp, indent, title, block);
192 for (i = 0; i <= block->nr; i++)
194 pr_indent(fp, indent+INDENT);
195 fprintf(fp, "%s->index[%d]=%d\n",
196 title, bShowNumbers ? i : -1, block->index[i]);
198 for (i = 0; i < block->nra; i++)
200 pr_indent(fp, indent+INDENT);
201 fprintf(fp, "%s->a[%d]=%d\n",
202 title, bShowNumbers ? i : -1, block->a[i]);
207 void pr_block(FILE *fp, int indent, const char *title, const t_block *block, gmx_bool bShowNumbers)
209 int i, start;
211 if (available(fp, block, indent, title))
213 indent = pr_block_title(fp, indent, title, block);
214 start = 0;
215 if (block->index[start] != 0)
217 fprintf(fp, "block->index[%d] should be 0\n", start);
219 else
221 for (i = 0; i < block->nr; i++)
223 int end = block->index[i+1];
224 pr_indent(fp, indent);
225 if (end <= start)
227 fprintf(fp, "%s[%d]={}\n", title, i);
229 else
231 fprintf(fp, "%s[%d]={%d..%d}\n",
232 title, bShowNumbers ? i : -1,
233 bShowNumbers ? start : -1, bShowNumbers ? end-1 : -1);
235 start = end;
241 void pr_blocka(FILE *fp, int indent, const char *title, const t_blocka *block, gmx_bool bShowNumbers)
243 int i, j, ok, size, start, end;
245 if (available(fp, block, indent, title))
247 indent = pr_blocka_title(fp, indent, title, block);
248 start = 0;
249 end = start;
250 if ((ok = (block->index[start] == 0)) == 0)
252 fprintf(fp, "block->index[%d] should be 0\n", start);
254 else
256 for (i = 0; i < block->nr; i++)
258 end = block->index[i+1];
259 size = pr_indent(fp, indent);
260 if (end <= start)
262 size += fprintf(fp, "%s[%d]={", title, i);
264 else
266 size += fprintf(fp, "%s[%d][%d..%d]={",
267 title, bShowNumbers ? i : -1,
268 bShowNumbers ? start : -1, bShowNumbers ? end-1 : -1);
270 for (j = start; j < end; j++)
272 if (j > start)
274 size += fprintf(fp, ", ");
276 if ((size) > (USE_WIDTH))
278 fprintf(fp, "\n");
279 size = pr_indent(fp, indent+INDENT);
281 size += fprintf(fp, "%d", block->a[j]);
283 fprintf(fp, "}\n");
284 start = end;
287 if ((end != block->nra) || (!ok))
289 pr_indent(fp, indent);
290 fprintf(fp, "tables inconsistent, dumping complete tables:\n");
291 low_pr_blocka(fp, indent, title, block, bShowNumbers);