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
45 #include <gmx_fatal.h>
48 #include "gromacs/selection/indexutil.h"
52 //! Alignment in bytes for all returned blocks.
56 * Describes a single block allocated from the memory pool.
58 typedef struct gmx_sel_mempool_block_t
60 //! Pointer to the start of the block (as returned to the user).
62 //! Size of the block, including padding required to align next block.
64 } gmx_sel_mempool_block_t;
67 * Describes a memory pool.
69 struct gmx_sel_mempool_t
71 //! Number of bytes currently allocated from the pool.
73 //! Number of bytes free in the pool, or 0 if \a buffer is NULL.
75 //! Memory area allocated for the pool, or NULL if not yet reserved.
77 //! Pointer to the first free byte (aligned at ::ALIGN_STEP) in \a buffer.
79 //! Number of blocks allocated from the pool.
81 //! Array describing the allocated blocks.
82 gmx_sel_mempool_block_t *blockstack;
83 //! Number of elements allocated for the \a blockstack array.
84 int blockstack_nalloc;
86 * Maximum number of bytes that have been reserved from the pool
93 _gmx_sel_mempool_create(gmx_sel_mempool_t **mpp)
95 gmx_sel_mempool_t *mp;
103 mp->blockstack = NULL;
104 mp->blockstack_nalloc = 0;
111 _gmx_sel_mempool_destroy(gmx_sel_mempool_t *mp)
117 for (i = 0; i < mp->nblocks; ++i)
119 sfree(mp->blockstack[i].ptr);
123 sfree(mp->blockstack);
128 _gmx_sel_mempool_alloc(gmx_sel_mempool_t *mp, void **ptrp, size_t size)
134 size_walign = ((size + ALIGN_STEP - 1) / ALIGN_STEP) * ALIGN_STEP;
137 if (mp->freesize < size)
139 gmx_bug("out of memory pool memory");
143 mp->freeptr += size_walign;
144 mp->freesize -= size_walign;
145 mp->currsize += size_walign;
152 gmx_mem("out of memory");
155 mp->currsize += size_walign;
156 if (mp->currsize > mp->maxsize)
158 mp->maxsize = mp->currsize;
162 if (mp->nblocks >= mp->blockstack_nalloc)
164 mp->blockstack_nalloc = mp->nblocks + 10;
165 srenew(mp->blockstack, mp->blockstack_nalloc);
167 mp->blockstack[mp->nblocks].ptr = ptr;
168 mp->blockstack[mp->nblocks].size = size_walign;
176 _gmx_sel_mempool_free(gmx_sel_mempool_t *mp, void *ptr)
184 assert(mp->nblocks > 0 && mp->blockstack[mp->nblocks - 1].ptr == ptr);
186 size = mp->blockstack[mp->nblocks].size;
187 mp->currsize -= size;
190 mp->freeptr = (char *)ptr;
191 mp->freesize += size;
200 _gmx_sel_mempool_reserve(gmx_sel_mempool_t *mp, size_t size)
202 assert(mp->nblocks == 0 && !mp->buffer);
207 mp->buffer = (char *)malloc(size);
210 gmx_mem("out of memory");
214 mp->freeptr = mp->buffer;
219 _gmx_sel_mempool_alloc_group(gmx_sel_mempool_t *mp, gmx_ana_index_t *g,
222 return _gmx_sel_mempool_alloc(mp, (void **)&g->index,
223 sizeof(*g->index)*isize);
227 _gmx_sel_mempool_free_group(gmx_sel_mempool_t *mp, gmx_ana_index_t *g)
229 _gmx_sel_mempool_free(mp, g->index);