Add gmx convert-trj
[gromacs.git] / src / gromacs / selection / mempool.cpp
blobd7c7cd3759f99e901e582a82c7c6bb8b423eac9b
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2010,2011,2012,2014,2017,2018, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
35 /*! \internal \file
36 * \brief
37 * Implements functions in mempool.h.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_selection
42 #include "gmxpre.h"
44 #include "mempool.h"
46 #include <cstdlib>
48 #include <new>
50 #include "gromacs/selection/indexutil.h"
51 #include "gromacs/utility/exceptions.h"
52 #include "gromacs/utility/gmxassert.h"
53 #include "gromacs/utility/smalloc.h"
55 //! Alignment in bytes for all returned blocks.
56 #define ALIGN_STEP 8
58 /*! \internal \brief
59 * Describes a single block allocated from the memory pool.
61 typedef struct gmx_sel_mempool_block_t
63 //! Pointer to the start of the block (as returned to the user).
64 void *ptr;
65 //! Size of the block, including padding required to align next block.
66 size_t size;
67 } gmx_sel_mempool_block_t;
69 /*! \internal \brief
70 * Describes a memory pool.
72 struct gmx_sel_mempool_t
74 //! Number of bytes currently allocated from the pool.
75 size_t currsize;
76 //! Number of bytes free in the pool, or 0 if \a buffer is NULL.
77 size_t freesize;
78 //! Memory area allocated for the pool, or NULL if not yet reserved.
79 char *buffer;
80 //! Pointer to the first free byte (aligned at ::ALIGN_STEP) in \a buffer.
81 char *freeptr;
82 //! Number of blocks allocated from the pool.
83 int nblocks;
84 //! Array describing the allocated blocks.
85 gmx_sel_mempool_block_t *blockstack;
86 //! Number of elements allocated for the \a blockstack array.
87 int blockstack_nalloc;
88 /*! \brief
89 * Maximum number of bytes that have been reserved from the pool
90 * simultaneously.
92 size_t maxsize;
95 gmx_sel_mempool_t *
96 _gmx_sel_mempool_create()
98 gmx_sel_mempool_t *mp;
100 snew(mp, 1);
101 mp->currsize = 0;
102 mp->freesize = 0;
103 mp->buffer = nullptr;
104 mp->freeptr = nullptr;
105 mp->nblocks = 0;
106 mp->blockstack = nullptr;
107 mp->blockstack_nalloc = 0;
108 mp->maxsize = 0;
109 return mp;
112 void
113 _gmx_sel_mempool_destroy(gmx_sel_mempool_t *mp)
115 if (!mp->buffer)
117 int i;
119 for (i = 0; i < mp->nblocks; ++i)
121 sfree(mp->blockstack[i].ptr);
124 sfree(mp->buffer);
125 sfree(mp->blockstack);
126 sfree(mp);
129 void *
130 _gmx_sel_mempool_alloc(gmx_sel_mempool_t *mp, size_t size)
132 void *ptr = nullptr;
133 size_t size_walign;
135 size_walign = ((size + ALIGN_STEP - 1) / ALIGN_STEP) * ALIGN_STEP;
136 if (mp->buffer)
138 if (mp->freesize < size)
140 GMX_THROW(gmx::InternalError("Out of memory pool memory"));
142 ptr = mp->freeptr;
143 mp->freeptr += size_walign;
144 mp->freesize -= size_walign;
145 mp->currsize += size_walign;
147 else
149 ptr = malloc(size);
150 if (!ptr)
152 throw std::bad_alloc();
154 mp->currsize += size_walign;
155 if (mp->currsize > mp->maxsize)
157 mp->maxsize = mp->currsize;
161 if (mp->nblocks >= mp->blockstack_nalloc)
163 mp->blockstack_nalloc = mp->nblocks + 10;
164 srenew(mp->blockstack, mp->blockstack_nalloc);
166 mp->blockstack[mp->nblocks].ptr = ptr;
167 mp->blockstack[mp->nblocks].size = size_walign;
168 mp->nblocks++;
170 return ptr;
173 void
174 _gmx_sel_mempool_free(gmx_sel_mempool_t *mp, void *ptr)
176 int size;
178 if (ptr == nullptr)
180 return;
182 GMX_RELEASE_ASSERT(mp->nblocks > 0 && mp->blockstack[mp->nblocks - 1].ptr == ptr,
183 "Invalid order of memory pool free calls");
184 mp->nblocks--;
185 size = mp->blockstack[mp->nblocks].size;
186 mp->currsize -= size;
187 if (mp->buffer)
189 mp->freeptr = static_cast<char *>(ptr);
190 mp->freesize += size;
192 else
194 sfree(ptr);
198 void
199 _gmx_sel_mempool_reserve(gmx_sel_mempool_t *mp, size_t size)
201 GMX_RELEASE_ASSERT(mp->nblocks == 0,
202 "Cannot reserve memory pool when there is something allocated");
203 GMX_RELEASE_ASSERT(!mp->buffer, "Cannot reserve memory pool twice");
204 if (size == 0)
206 size = mp->maxsize;
208 mp->buffer = static_cast<char *>(malloc(size));
209 if (!mp->buffer)
211 throw std::bad_alloc();
213 mp->freesize = size;
214 mp->freeptr = mp->buffer;
217 void
218 _gmx_sel_mempool_alloc_group(gmx_sel_mempool_t *mp, gmx_ana_index_t *g,
219 int isize)
221 void *ptr = _gmx_sel_mempool_alloc(mp, sizeof(*g->index)*isize);
222 g->index = static_cast<int *>(ptr);
225 void
226 _gmx_sel_mempool_free_group(gmx_sel_mempool_t *mp, gmx_ana_index_t *g)
228 _gmx_sel_mempool_free(mp, g->index);
229 g->index = nullptr;