2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2010,2011,2012,2014,2017 by the GROMACS development team.
5 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
38 * Implements functions in mempool.h.
40 * \author Teemu Murtola <teemu.murtola@gmail.com>
41 * \ingroup module_selection
51 #include "gromacs/selection/indexutil.h"
52 #include "gromacs/utility/exceptions.h"
53 #include "gromacs/utility/gmxassert.h"
54 #include "gromacs/utility/smalloc.h"
56 //! Alignment in bytes for all returned blocks.
60 * Describes a single block allocated from the memory pool.
62 typedef struct gmx_sel_mempool_block_t
64 //! Pointer to the start of the block (as returned to the user).
66 //! Size of the block, including padding required to align next block.
68 } gmx_sel_mempool_block_t
;
71 * Describes a memory pool.
73 struct gmx_sel_mempool_t
75 //! Number of bytes currently allocated from the pool.
77 //! Number of bytes free in the pool, or 0 if \a buffer is NULL.
79 //! Memory area allocated for the pool, or NULL if not yet reserved.
81 //! Pointer to the first free byte (aligned at ::ALIGN_STEP) in \a buffer.
83 //! Number of blocks allocated from the pool.
85 //! Array describing the allocated blocks.
86 gmx_sel_mempool_block_t
* blockstack
;
87 //! Number of elements allocated for the \a blockstack array.
88 int blockstack_nalloc
;
90 * Maximum number of bytes that have been reserved from the pool
96 gmx_sel_mempool_t
* _gmx_sel_mempool_create()
98 gmx_sel_mempool_t
* mp
;
103 mp
->buffer
= nullptr;
104 mp
->freeptr
= nullptr;
106 mp
->blockstack
= nullptr;
107 mp
->blockstack_nalloc
= 0;
112 void _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
);
128 void* _gmx_sel_mempool_alloc(gmx_sel_mempool_t
* mp
, size_t size
)
133 size_walign
= ((size
+ ALIGN_STEP
- 1) / ALIGN_STEP
) * ALIGN_STEP
;
136 if (mp
->freesize
< size
)
138 GMX_THROW(gmx::InternalError("Out of memory pool memory"));
141 mp
->freeptr
+= size_walign
;
142 mp
->freesize
-= size_walign
;
143 mp
->currsize
+= size_walign
;
150 throw std::bad_alloc();
152 mp
->currsize
+= size_walign
;
153 if (mp
->currsize
> mp
->maxsize
)
155 mp
->maxsize
= mp
->currsize
;
159 if (mp
->nblocks
>= mp
->blockstack_nalloc
)
161 mp
->blockstack_nalloc
= mp
->nblocks
+ 10;
162 srenew(mp
->blockstack
, mp
->blockstack_nalloc
);
164 mp
->blockstack
[mp
->nblocks
].ptr
= ptr
;
165 mp
->blockstack
[mp
->nblocks
].size
= size_walign
;
171 void _gmx_sel_mempool_free(gmx_sel_mempool_t
* mp
, void* ptr
)
179 GMX_RELEASE_ASSERT(mp
->nblocks
> 0 && mp
->blockstack
[mp
->nblocks
- 1].ptr
== ptr
,
180 "Invalid order of memory pool free calls");
182 size
= mp
->blockstack
[mp
->nblocks
].size
;
183 mp
->currsize
-= size
;
186 mp
->freeptr
= static_cast<char*>(ptr
);
187 mp
->freesize
+= size
;
195 void _gmx_sel_mempool_reserve(gmx_sel_mempool_t
* mp
, size_t size
)
197 GMX_RELEASE_ASSERT(mp
->nblocks
== 0,
198 "Cannot reserve memory pool when there is something allocated");
199 GMX_RELEASE_ASSERT(!mp
->buffer
, "Cannot reserve memory pool twice");
204 mp
->buffer
= static_cast<char*>(malloc(size
));
207 throw std::bad_alloc();
210 mp
->freeptr
= mp
->buffer
;
213 void _gmx_sel_mempool_alloc_group(gmx_sel_mempool_t
* mp
, gmx_ana_index_t
* g
, int isize
)
215 void* ptr
= _gmx_sel_mempool_alloc(mp
, sizeof(*g
->index
) * isize
);
216 g
->index
= static_cast<int*>(ptr
);
219 void _gmx_sel_mempool_free_group(gmx_sel_mempool_t
* mp
, gmx_ana_index_t
* g
)
221 _gmx_sel_mempool_free(mp
, g
->index
);