2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2018,2019,2020, by the GROMACS development team.
7 * Copyright (c) 2021, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
40 * C-style memory allocation routines for \Gromacs.
42 * This header provides macros snew(), srenew(), smalloc(), and sfree() for
43 * C-style memory management. Additionally, snew_aligned() and sfree_aligned() are
44 * provided for managing memory with a specified byte alignment.
46 * If an allocation fails, the program is halted by calling gmx_fatal(), which
47 * outputs source file and line number and the name of the variable involved.
48 * This frees calling code from the trouble of checking the result of the
49 * allocations everywhere. It also provides a location for centrally logging
50 * memory allocations for diagnosing memory usage (currently can only enabled
51 * by changing the source code). Additionally, sfree() works also with a
52 * `NULL` parameter, which standard free() does not.
54 * The macros forward the calls to functions save_malloc(), save_calloc(),
55 * save_realloc(), save_free(), save_calloc_aligned(), and save_free_aligned().
56 * There are a few low-level locations in \Gromacs that call these directly,
57 * but generally the macros should be used.
58 * save_malloc_aligned() exists for this purpose, although there is no macro to
62 * \ingroup module_utility
64 #ifndef GMX_UTILITY_SMALLOC_H
65 #define GMX_UTILITY_SMALLOC_H
70 * \Gromacs wrapper for malloc().
72 * \param[in] name Variable name identifying the allocation.
73 * \param[in] file Source code file where the allocation originates from.
74 * \param[in] line Source code line where the allocation originates from.
75 * \param[in] size Number of bytes to allocate.
76 * \returns Pointer to the allocated space.
78 * This should generally be called through smalloc(), not directly.
80 void* save_malloc(const char* name, const char* file, int line, size_t size);
82 * \Gromacs wrapper for calloc().
84 * \param[in] name Variable name identifying the allocation.
85 * \param[in] file Source code file where the allocation originates from.
86 * \param[in] line Source code line where the allocation originates from.
87 * \param[in] nelem Number of elements to allocate.
88 * \param[in] elsize Number of bytes per element.
89 * \returns Pointer to the allocated space.
91 * This should generally be called through snew(), not directly.
93 void* save_calloc(const char* name, const char* file, int line, size_t nelem, size_t elsize);
95 * \Gromacs wrapper for realloc().
97 * \param[in] name Variable name identifying the allocation.
98 * \param[in] file Source code file where the allocation originates from.
99 * \param[in] line Source code line where the allocation originates from.
100 * \param[in] ptr Pointer to the previously allocated memory (can be NULL).
101 * \param[in] nelem Number of elements to allocate.
102 * \param[in] elsize Number of bytes per element.
103 * \returns Pointer to the allocated space.
105 * As with realloc(), if \p ptr is NULL, memory is allocated as if malloc() was
107 * This should generally be called through srenew(), not directly.
109 * Note that the allocated memory is not initialized to zero.
111 void* save_realloc(const char* name, const char* file, int line, void* ptr, size_t nelem, size_t elsize);
113 * \Gromacs wrapper for free().
115 * \param[in] name Variable name identifying the deallocation.
116 * \param[in] file Source code file where the deallocation originates from.
117 * \param[in] line Source code line where the deallocation originates from.
118 * \param[in] ptr Pointer to the allocated memory (can be NULL).
120 * If \p ptr is NULL, does nothing.
121 * This should generally be called through sfree(), not directly.
124 void save_free(const char* name, const char* file, int line, void* ptr);
127 * \Gromacs wrapper for allocating aligned memory.
129 * \param[in] name Variable name identifying the allocation.
130 * \param[in] file Source code file where the allocation originates from.
131 * \param[in] line Source code line where the allocation originates from.
132 * \param[in] nelem Number of elements to allocate.
133 * \param[in] elsize Number of bytes per element.
134 * \param[in] alignment Requested alignment in bytes.
135 * \returns Pointer to the allocated space, aligned at `alignment`-byte
138 * There is no macro that invokes this function.
140 * The returned pointer should only be freed with a call to save_free_aligned().
142 void* save_malloc_aligned(const char* name, const char* file, int line, size_t nelem, size_t elsize, size_t alignment);
144 * \Gromacs wrapper for allocating zero-initialized aligned memory.
146 * \param[in] name Variable name identifying the allocation.
147 * \param[in] file Source code file where the allocation originates from.
148 * \param[in] line Source code line where the allocation originates from.
149 * \param[in] nelem Number of elements to allocate.
150 * \param[in] elsize Number of bytes per element.
151 * \param[in] alignment Requested alignment in bytes.
152 * \returns Pointer to the allocated space, aligned at `alignment`-byte
155 * This should generally be called through snew_aligned(), not directly.
157 * The returned pointer should only be freed with a call to save_free_aligned().
159 void* save_calloc_aligned(const char* name, const char* file, int line, size_t nelem, size_t elsize, size_t alignment);
161 * \Gromacs wrapper for freeing aligned memory.
163 * \param[in] name Variable name identifying the deallocation.
164 * \param[in] file Source code file where the deallocation originates from.
165 * \param[in] line Source code line where the deallocation originates from.
166 * \param[in] ptr Pointer to the allocated memory (can be NULL).
168 * If \p ptr is NULL, does nothing.
169 * \p ptr should have been allocated with save_malloc_aligned() or
170 * save_calloc_aligned().
171 * This should generally be called through sfree_aligned(), not directly.
174 void save_free_aligned(const char* name, const char* file, int line, void* ptr);
176 #include <type_traits>
178 /*! \cond internal */
179 /*! \name Implementation templates for C++ memory allocation macros
181 * These templates are used to implement the snew() etc. macros for C++, where
182 * an explicit cast is needed from `void *` (the return value of the allocation
183 * wrapper functions) to the type of \p ptr.
185 * Having these as `static` avoid some obscure bugs if several files define
186 * distinct data structures with identical names and allocate memory for them
187 * using snew(). By the C++ standard, such declarations cause undefined
188 * behavior, but can be difficult to spot in the existing C code.
189 * Without the `static` (and if the compiler does not inline the calls), the
190 * linker cannot that data structures with identical names are actually
191 * different and links calls to these template functions incorrectly, which can
192 * result in allocation of an incorrect amount of memory if the element size is
193 * computed within the function.
195 * The size cannot be passed as a parameter to the function either, since that
196 * provokes warnings from cppcheck for some invocations, where a complex
197 * expression is passed as \p ptr.
200 /** C++ helper for snew(). */
202 static inline void gmx_snew_impl(const char* name, const char* file, int line, T*& ptr, size_t nelem)
204 // TODO: Use std::is_pod_v when CUDA 11 is a requirement.
205 static_assert(std::is_pod<T>::value, "snew() called on C++ type");
206 // NOLINTNEXTLINE bugprone-sizeof-expression
207 ptr = static_cast<T*>(save_calloc(name, file, line, nelem, sizeof(T)));
209 /** C++ helper for srenew(). */
211 static inline void gmx_srenew_impl(const char* name, const char* file, int line, T*& ptr, size_t nelem)
213 // TODO: Use std::is_pod_v when CUDA 11 is a requirement.
214 static_assert(std::is_pod<T>::value, "srenew() called on C++ type");
215 // NOLINTNEXTLINE bugprone-sizeof-expression
216 ptr = static_cast<T*>(save_realloc(name, file, line, ptr, nelem, sizeof(T)));
218 /** C++ helper for smalloc(). */
220 static inline void gmx_smalloc_impl(const char* name, const char* file, int line, T*& ptr, size_t size)
222 // TODO: Use std::is_pod_v when CUDA 11 is a requirement.
223 static_assert(std::is_pod<T>::value, "smalloc() called on C++ type");
224 ptr = static_cast<T*>(save_malloc(name, file, line, size));
226 /** C++ helper for snew_aligned(). */
229 gmx_snew_aligned_impl(const char* name, const char* file, int line, T*& ptr, size_t nelem, size_t alignment)
231 // TODO: Use std::is_pod_v when CUDA 11 is a requirement.
232 static_assert(std::is_pod<T>::value, "snew_aligned() called on C++ type");
233 ptr = static_cast<T*>(save_calloc_aligned(name, file, line, nelem, sizeof(T), alignment));
235 /** C++ helper for sfree(). */
237 static inline void gmx_sfree_impl(const char* name, const char* file, int line, T* ptr)
239 // TODO: Use std::is_pod_v and std::is_void_v when CUDA 11 is a requirement.
240 static_assert(std::is_pod<T>::value || std::is_void<T>::value, "sfree() called on C++ type");
241 save_free(name, file, line, ptr);
243 /** C++ helper for sfree_aligned(). */
245 static inline void gmx_sfree_aligned_impl(const char* name, const char* file, int line, T* ptr)
247 // TODO: Use std::is_pod_v and std::is_void_v when CUDA 11 is a requirement.
248 static_assert(std::is_pod<T>::value || std::is_void<T>::value,
249 "sfree_aligned() called on C++ type");
250 save_free_aligned(name, file, line, ptr);
257 * Allocates memory for a given number of elements.
259 * \param[out] ptr Pointer to allocate.
260 * \param[in] nelem Number of elements to allocate.
262 * Allocates memory for \p nelem elements of type \p *ptr and sets this to
263 * \p ptr. The allocated memory is initialized to zeros.
269 * Reallocates memory for a given number of elements.
271 * \param[in,out] ptr Pointer to allocate/reallocate.
272 * \param[in] nelem Number of elements to allocate.
274 * (Re)allocates memory for \p ptr such that it can hold \p nelem elements of
275 * type \p *ptr, and sets the new pointer to \p ptr.
276 * If \p ptr is `NULL`, memory is allocated as if it was new.
277 * If \p nelem is zero, \p ptr is freed (if not `NULL`).
278 * Note that the allocated memory is not initialized, unlike with snew().
284 * Allocates memory for a given number of bytes.
286 * \param[out] ptr Pointer to allocate.
287 * \param[in] size Number of bytes to allocate.
289 * Allocates memory for \p size bytes and sets this to \p ptr.
290 * The allocated memory is initialized to zero.
294 /*! \def snew_aligned
296 * Allocates aligned memory for a given number of elements.
298 * \param[out] ptr Pointer to allocate.
299 * \param[in] nelem Number of elements to allocate.
300 * \param[in] alignment Requested alignment in bytes.
302 * Allocates memory for \p nelem elements of type \p *ptr and sets this to
303 * \p ptr. The returned pointer is `alignment`-byte aligned.
304 * The allocated memory is initialized to zeros.
306 * The returned pointer should only be freed with sfree_aligned().
312 * Frees memory referenced by \p ptr.
314 * \p ptr is allowed to be NULL, in which case nothing is done.
318 /*! \def sfree_aligned
320 * Frees aligned memory referenced by \p ptr.
322 * This must only be called with a pointer obtained through snew_aligned().
323 * \p ptr is allowed to be NULL, in which case nothing is done.
328 /* C++ implementation */
329 #define snew(ptr, nelem) gmx_snew_impl(#ptr, __FILE__, __LINE__, (ptr), (nelem))
330 #define srenew(ptr, nelem) gmx_srenew_impl(#ptr, __FILE__, __LINE__, (ptr), (nelem))
331 #define smalloc(ptr, size) gmx_smalloc_impl(#ptr, __FILE__, __LINE__, (ptr), (size))
332 #define snew_aligned(ptr, nelem, alignment) \
333 gmx_snew_aligned_impl(#ptr, __FILE__, __LINE__, (ptr), (nelem), alignment)
334 #define sfree(ptr) gmx_sfree_impl(#ptr, __FILE__, __LINE__, (ptr))
335 #define sfree_aligned(ptr) gmx_sfree_aligned_impl(#ptr, __FILE__, __LINE__, (ptr))
338 * Over allocation factor for memory allocations.
340 * Memory (re)allocation can be VERY slow, especially with some
341 * MPI libraries that replace the standard malloc and realloc calls.
342 * To avoid slow memory allocation we use over_alloc to set the memory
343 * allocation size for large data blocks. Since this scales the size
344 * with a factor, we use log(n) realloc calls instead of n.
345 * This can reduce allocation times from minutes to seconds.
347 * This factor leads to 4 realloc calls to double the array size.
349 constexpr float OVER_ALLOC_FAC = 1.19F;
352 * Turns over allocation for variable size atoms/cg/top arrays on or off,
356 * This is mdrun-specific, so it might be better to put this and
357 * over_alloc_dd() much higher up.
359 void set_over_alloc_dd(bool set);
362 * Returns new allocation count for domain decomposition allocations.
364 * Returns n when domain decomposition over allocation is off.
365 * Returns OVER_ALLOC_FAC*n + 100 when over allocation in on.
366 * This is to avoid frequent reallocation during domain decomposition in mdrun.
368 int over_alloc_dd(int n);
370 /** Over allocation for small data types: int, real etc. */
372 constexpr T over_alloc_small(T n)
374 return OVER_ALLOC_FAC * n + 8000;
377 /** Over allocation for large data types: complex structs */
379 constexpr T over_alloc_large(T n)
381 return OVER_ALLOC_FAC * n + 1000;