Merge release-4-6 (commit 'Ic142a690')
[alexxy/gromacs.git] / src / gromacs / selection / mempool.cpp
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
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.
13
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.
18  *
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.
25  *
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.
28  *
29  * For more info, check our website at http://www.gromacs.org
30  */
31 /*! \internal \file
32  * \brief
33  * Implements functions in mempool.h.
34  *
35  * \author Teemu Murtola <teemu.murtola@cbr.su.se>
36  * \ingroup module_selection
37  */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <stdlib.h>
43
44 #include <new>
45
46 #include "smalloc.h"
47
48 #include "gromacs/selection/indexutil.h"
49 #include "gromacs/utility/exceptions.h"
50 #include "gromacs/utility/gmxassert.h"
51
52 #include "mempool.h"
53
54 //! Alignment in bytes for all returned blocks.
55 #define ALIGN_STEP 8
56
57 /*! \internal \brief
58  * Describes a single block allocated from the memory pool.
59  */
60 typedef struct gmx_sel_mempool_block_t
61 {
62     //! Pointer to the start of the block (as returned to the user).
63     void                       *ptr;
64     //! Size of the block, including padding required to align next block.
65     size_t                      size;
66 } gmx_sel_mempool_block_t;
67
68 /*! \internal \brief
69  * Describes a memory pool.
70  */
71 struct gmx_sel_mempool_t
72 {
73     //! Number of bytes currently allocated from the pool.
74     size_t                      currsize;
75     //! Number of bytes free in the pool, or 0 if \a buffer is NULL.
76     size_t                      freesize;
77     //! Memory area allocated for the pool, or NULL if not yet reserved.
78     char                       *buffer;
79     //! Pointer to the first free byte (aligned at ::ALIGN_STEP) in \a buffer.
80     char                       *freeptr;
81     //! Number of blocks allocated from the pool.
82     int                         nblocks;
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;
87     /*! \brief
88      * Maximum number of bytes that have been reserved from the pool
89      * simultaneously.
90      */
91     size_t                      maxsize;
92 };
93
94 gmx_sel_mempool_t *
95 _gmx_sel_mempool_create()
96 {
97     gmx_sel_mempool_t *mp;
98
99     snew(mp, 1);
100     mp->currsize          = 0;
101     mp->freesize          = 0;
102     mp->buffer            = NULL;
103     mp->freeptr           = NULL;
104     mp->nblocks           = 0;
105     mp->blockstack        = NULL;
106     mp->blockstack_nalloc = 0;
107     mp->maxsize           = 0;
108     return mp;
109 }
110
111 void
112 _gmx_sel_mempool_destroy(gmx_sel_mempool_t *mp)
113 {
114     if (!mp->buffer)
115     {
116         int  i;
117
118         for (i = 0; i < mp->nblocks; ++i)
119         {
120             sfree(mp->blockstack[i].ptr);
121         }
122     }
123     sfree(mp->buffer);
124     sfree(mp->blockstack);
125     sfree(mp);
126 }
127
128 void *
129 _gmx_sel_mempool_alloc(gmx_sel_mempool_t *mp, size_t size)
130 {
131     void   *ptr = NULL;
132     size_t  size_walign;
133
134     size_walign = ((size + ALIGN_STEP - 1) / ALIGN_STEP) * ALIGN_STEP;
135     if (mp->buffer)
136     {
137         if (mp->freesize < size)
138         {
139             GMX_THROW(gmx::InternalError("Out of memory pool memory"));
140         }
141         ptr = mp->freeptr;
142         mp->freeptr  += size_walign;
143         mp->freesize -= size_walign;
144         mp->currsize += size_walign;
145     }
146     else
147     {
148         ptr = malloc(size);
149         if (!ptr)
150         {
151             throw std::bad_alloc();
152         }
153         mp->currsize += size_walign;
154         if (mp->currsize > mp->maxsize)
155         {
156             mp->maxsize = mp->currsize;
157         }
158     }
159
160     if (mp->nblocks >= mp->blockstack_nalloc)
161     {
162         mp->blockstack_nalloc = mp->nblocks + 10;
163         srenew(mp->blockstack, mp->blockstack_nalloc);
164     }
165     mp->blockstack[mp->nblocks].ptr  = ptr;
166     mp->blockstack[mp->nblocks].size = size_walign;
167     mp->nblocks++;
168
169     return ptr;
170 }
171
172 void
173 _gmx_sel_mempool_free(gmx_sel_mempool_t *mp, void *ptr)
174 {
175     int size;
176
177     if (ptr == NULL)
178     {
179         return;
180     }
181     GMX_RELEASE_ASSERT(mp->nblocks > 0 && mp->blockstack[mp->nblocks - 1].ptr == ptr,
182                        "Invalid order of memory pool free calls");
183     mp->nblocks--;
184     size = mp->blockstack[mp->nblocks].size;
185     mp->currsize -= size;
186     if (mp->buffer)
187     {
188         mp->freeptr = (char *)ptr;
189         mp->freesize += size;
190     }
191     else
192     {
193         sfree(ptr);
194     }
195 }
196
197 void
198 _gmx_sel_mempool_reserve(gmx_sel_mempool_t *mp, size_t size)
199 {
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");
203     if (size == 0)
204     {
205         size = mp->maxsize;
206     }
207     mp->buffer = (char *)malloc(size);
208     if (!mp->buffer)
209     {
210         throw std::bad_alloc();
211     }
212     mp->freesize = size;
213     mp->freeptr  = mp->buffer;
214 }
215
216 void
217 _gmx_sel_mempool_alloc_group(gmx_sel_mempool_t *mp, gmx_ana_index_t *g,
218                              int isize)
219 {
220     void *ptr = _gmx_sel_mempool_alloc(mp, sizeof(*g->index)*isize);
221     g->index = static_cast<int *>(ptr);
222 }
223
224 void
225 _gmx_sel_mempool_free_group(gmx_sel_mempool_t *mp, gmx_ana_index_t *g)
226 {
227     _gmx_sel_mempool_free(mp, g->index);
228     g->index = NULL;
229 }