3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
9 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11 * Copyright (c) 2001-2009, The GROMACS development team,
12 * check out http://www.gromacs.org for more information.
14 * This program is free software; you can redistribute it and/or
15 * modify it under the terms of the GNU General Public License
16 * as published by the Free Software Foundation; either version 2
17 * of the License, or (at your option) any later version.
19 * If you want to redistribute modifications, please consider that
20 * scientific software is very special. Version control is crucial -
21 * bugs must be traceable. We will be happy to consider code for
22 * inclusion in the official distribution, but derived work must not
23 * be called official GROMACS. Details are found in the README & COPYING
24 * files - if they are missing, get the official version at www.gromacs.org.
26 * To help us fund GROMACS development, we humbly ask that you cite
27 * the papers on the package - you can find them in the top README file.
29 * For more info, check our website at http://www.gromacs.org
33 * Implements functions in mempool.h.
35 * \author Teemu Murtola <teemu.murtola@cbr.su.se>
36 * \ingroup module_selection
48 #include "gromacs/fatalerror/exceptions.h"
49 #include "gromacs/fatalerror/gmxassert.h"
50 #include "gromacs/selection/indexutil.h"
54 //! Alignment in bytes for all returned blocks.
58 * Describes a single block allocated from the memory pool.
60 typedef struct gmx_sel_mempool_block_t
62 //! Pointer to the start of the block (as returned to the user).
64 //! Size of the block, including padding required to align next block.
66 } gmx_sel_mempool_block_t;
69 * Describes a memory pool.
71 struct gmx_sel_mempool_t
73 //! Number of bytes currently allocated from the pool.
75 //! Number of bytes free in the pool, or 0 if \a buffer is NULL.
77 //! Memory area allocated for the pool, or NULL if not yet reserved.
79 //! Pointer to the first free byte (aligned at ::ALIGN_STEP) in \a buffer.
81 //! Number of blocks allocated from the pool.
83 //! Array describing the allocated blocks.
84 gmx_sel_mempool_block_t *blockstack;
85 //! Number of elements allocated for the \a blockstack array.
86 int blockstack_nalloc;
88 * Maximum number of bytes that have been reserved from the pool
95 _gmx_sel_mempool_create()
97 gmx_sel_mempool_t *mp;
105 mp->blockstack = NULL;
106 mp->blockstack_nalloc = 0;
112 _gmx_sel_mempool_destroy(gmx_sel_mempool_t *mp)
118 for (i = 0; i < mp->nblocks; ++i)
120 sfree(mp->blockstack[i].ptr);
124 sfree(mp->blockstack);
129 _gmx_sel_mempool_alloc(gmx_sel_mempool_t *mp, size_t size)
134 size_walign = ((size + ALIGN_STEP - 1) / ALIGN_STEP) * ALIGN_STEP;
137 if (mp->freesize < size)
139 GMX_THROW(gmx::InternalError("Out of memory pool memory"));
142 mp->freeptr += size_walign;
143 mp->freesize -= size_walign;
144 mp->currsize += size_walign;
151 throw std::bad_alloc();
153 mp->currsize += size_walign;
154 if (mp->currsize > mp->maxsize)
156 mp->maxsize = mp->currsize;
160 if (mp->nblocks >= mp->blockstack_nalloc)
162 mp->blockstack_nalloc = mp->nblocks + 10;
163 srenew(mp->blockstack, mp->blockstack_nalloc);
165 mp->blockstack[mp->nblocks].ptr = ptr;
166 mp->blockstack[mp->nblocks].size = size_walign;
173 _gmx_sel_mempool_free(gmx_sel_mempool_t *mp, void *ptr)
181 GMX_RELEASE_ASSERT(mp->nblocks > 0 && mp->blockstack[mp->nblocks - 1].ptr == ptr,
182 "Invalid order of memory pool free calls");
184 size = mp->blockstack[mp->nblocks].size;
185 mp->currsize -= size;
188 mp->freeptr = (char *)ptr;
189 mp->freesize += size;
198 _gmx_sel_mempool_reserve(gmx_sel_mempool_t *mp, size_t size)
200 GMX_RELEASE_ASSERT(mp->nblocks == 0,
201 "Cannot reserve memory pool when there is something allocated");
202 GMX_RELEASE_ASSERT(!mp->buffer, "Cannot reserve memory pool twice");
207 mp->buffer = (char *)malloc(size);
210 throw std::bad_alloc();
213 mp->freeptr = mp->buffer;
217 _gmx_sel_mempool_alloc_group(gmx_sel_mempool_t *mp, gmx_ana_index_t *g,
220 void *ptr = _gmx_sel_mempool_alloc(mp, sizeof(*g->index)*isize);
221 g->index = static_cast<int *>(ptr);
225 _gmx_sel_mempool_free_group(gmx_sel_mempool_t *mp, gmx_ana_index_t *g)
227 _gmx_sel_mempool_free(mp, g->index);