More fixes from clang-tidy-7
[gromacs.git] / src / gromacs / utility / smalloc.h
blob906da1c66c9220f818a12629b15a6687f73cf0ac
1 /*
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, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 /*! \file
38 * \brief
39 * C-style memory allocation routines for \Gromacs.
41 * This header provides macros snew(), srenew(), smalloc(), and sfree() for
42 * C-style memory management. Additionally, snew_aligned() and sfree_aligned() are
43 * provided for managing memory with a specified byte alignment.
45 * If an allocation fails, the program is halted by calling gmx_fatal(), which
46 * outputs source file and line number and the name of the variable involved.
47 * This frees calling code from the trouble of checking the result of the
48 * allocations everywhere. It also provides a location for centrally logging
49 * memory allocations for diagnosing memory usage (currently can only enabled
50 * by changing the source code). Additionally, sfree() works also with a
51 * `NULL` parameter, which standard free() does not.
53 * The macros forward the calls to functions save_malloc(), save_calloc(),
54 * save_realloc(), save_free(), save_calloc_aligned(), and save_free_aligned().
55 * There are a few low-level locations in \Gromacs that call these directly,
56 * but generally the macros should be used.
57 * save_malloc_aligned() exists for this purpose, although there is no macro to
58 * invoke it.
60 * \inpublicapi
61 * \ingroup module_utility
63 #ifndef GMX_UTILITY_SMALLOC_H
64 #define GMX_UTILITY_SMALLOC_H
66 #include <stddef.h>
68 #include "gromacs/utility/basedefinitions.h"
70 /*! \brief
71 * \Gromacs wrapper for malloc().
73 * \param[in] name Variable name identifying the allocation.
74 * \param[in] file Source code file where the allocation originates from.
75 * \param[in] line Source code line where the allocation originates from.
76 * \param[in] size Number of bytes to allocate.
77 * \returns Pointer to the allocated space.
79 * This should generally be called through smalloc(), not directly.
81 void *save_malloc(const char *name, const char *file, int line, size_t size);
82 /*! \brief
83 * \Gromacs wrapper for calloc().
85 * \param[in] name Variable name identifying the allocation.
86 * \param[in] file Source code file where the allocation originates from.
87 * \param[in] line Source code line where the allocation originates from.
88 * \param[in] nelem Number of elements to allocate.
89 * \param[in] elsize Number of bytes per element.
90 * \returns Pointer to the allocated space.
92 * This should generally be called through snew(), not directly.
94 void *save_calloc(const char *name, const char *file, int line,
95 size_t nelem, size_t elsize);
96 /*! \brief
97 * \Gromacs wrapper for realloc().
99 * \param[in] name Variable name identifying the allocation.
100 * \param[in] file Source code file where the allocation originates from.
101 * \param[in] line Source code line where the allocation originates from.
102 * \param[in] ptr Pointer to the previously allocated memory (can be NULL).
103 * \param[in] nelem Number of elements to allocate.
104 * \param[in] elsize Number of bytes per element.
105 * \returns Pointer to the allocated space.
107 * As with realloc(), if \p ptr is NULL, memory is allocated as if malloc() was
108 * called.
109 * This should generally be called through srenew(), not directly.
111 * Note that the allocated memory is not initialized to zero.
113 void *save_realloc(const char *name, const char *file, int line,
114 void *ptr, size_t nelem, size_t elsize);
115 /*! \brief
116 * \Gromacs wrapper for free().
118 * \param[in] name Variable name identifying the deallocation.
119 * \param[in] file Source code file where the deallocation originates from.
120 * \param[in] line Source code line where the deallocation originates from.
121 * \param[in] ptr Pointer to the allocated memory (can be NULL).
123 * If \p ptr is NULL, does nothing.
124 * This should generally be called through sfree(), not directly.
125 * This never fails.
127 void save_free(const char *name, const char *file, int line, void *ptr);
129 /*! \brief
130 * \Gromacs wrapper for allocating aligned memory.
132 * \param[in] name Variable name identifying the allocation.
133 * \param[in] file Source code file where the allocation originates from.
134 * \param[in] line Source code line where the allocation originates from.
135 * \param[in] nelem Number of elements to allocate.
136 * \param[in] elsize Number of bytes per element.
137 * \param[in] alignment Requested alignment in bytes.
138 * \returns Pointer to the allocated space, aligned at `alignment`-byte
139 * boundary.
141 * There is no macro that invokes this function.
143 * The returned pointer should only be freed with a call to save_free_aligned().
145 void *save_malloc_aligned(const char *name, const char *file, int line,
146 size_t nelem, size_t elsize, size_t alignment);
147 /*! \brief
148 * \Gromacs wrapper for allocating zero-initialized aligned memory.
150 * \param[in] name Variable name identifying the allocation.
151 * \param[in] file Source code file where the allocation originates from.
152 * \param[in] line Source code line where the allocation originates from.
153 * \param[in] nelem Number of elements to allocate.
154 * \param[in] elsize Number of bytes per element.
155 * \param[in] alignment Requested alignment in bytes.
156 * \returns Pointer to the allocated space, aligned at `alignment`-byte
157 * boundary.
159 * This should generally be called through snew_aligned(), not directly.
161 * The returned pointer should only be freed with a call to save_free_aligned().
163 void *save_calloc_aligned(const char *name, const char *file, int line,
164 size_t nelem, size_t elsize, size_t alignment);
165 /*! \brief
166 * \Gromacs wrapper for freeing aligned memory.
168 * \param[in] name Variable name identifying the deallocation.
169 * \param[in] file Source code file where the deallocation originates from.
170 * \param[in] line Source code line where the deallocation originates from.
171 * \param[in] ptr Pointer to the allocated memory (can be NULL).
173 * If \p ptr is NULL, does nothing.
174 * \p ptr should have been allocated with save_malloc_aligned() or
175 * save_calloc_aligned().
176 * This should generally be called through sfree_aligned(), not directly.
177 * This never fails.
179 void save_free_aligned(const char *name, const char *file, int line, void *ptr);
181 #include <type_traits>
183 /*! \cond internal */
184 /*! \name Implementation templates for C++ memory allocation macros
186 * These templates are used to implement the snew() etc. macros for C++, where
187 * an explicit cast is needed from `void *` (the return value of the allocation
188 * wrapper functions) to the type of \p ptr.
190 * Having these as `static` avoid some obscure bugs if several files define
191 * distinct data structures with identical names and allocate memory for them
192 * using snew(). By the C++ standard, such declarations cause undefined
193 * behavior, but can be difficult to spot in the existing C code.
194 * Without the `static` (and if the compiler does not inline the calls), the
195 * linker cannot that data structures with identical names are actually
196 * different and links calls to these template functions incorrectly, which can
197 * result in allocation of an incorrect amount of memory if the element size is
198 * computed within the function.
200 * The size cannot be passed as a parameter to the function either, since that
201 * provokes warnings from cppcheck for some invocations, where a complex
202 * expression is passed as \p ptr.
204 /*! \{ */
205 /** C++ helper for snew(). */
206 template <typename T> static inline
207 void gmx_snew_impl(const char *name, const char *file, int line,
208 T * &ptr, size_t nelem)
210 static_assert(std::is_pod<T>::value, "snew() called on C++ type");
211 ptr = static_cast<T*>(save_calloc(name, file, line, nelem, sizeof(T)));
213 /** C++ helper for srenew(). */
214 template <typename T> static inline
215 void gmx_srenew_impl(const char *name, const char *file, int line,
216 T * &ptr, size_t nelem)
218 static_assert(std::is_pod<T>::value, "srenew() called on C++ type");
219 ptr = static_cast<T*>(save_realloc(name, file, line, ptr, nelem, sizeof(T)));
221 /** C++ helper for smalloc(). */
222 template <typename T> static inline
223 void gmx_smalloc_impl(const char *name, const char *file, int line,
224 T * &ptr, size_t size)
226 static_assert(std::is_pod<T>::value, "smalloc() called on C++ type");
227 ptr = static_cast<T*>(save_malloc(name, file, line, size));
229 /** C++ helper for snew_aligned(). */
230 template <typename T> static inline
231 void gmx_snew_aligned_impl(const char *name, const char *file, int line,
232 T * &ptr, size_t nelem, size_t alignment)
234 static_assert(std::is_pod<T>::value, "snew_aligned() called on C++ type");
235 ptr = static_cast<T*>(save_calloc_aligned(name, file, line, nelem, sizeof(T), alignment));
237 /** C++ helper for sfree(). */
238 template <typename T> static inline
239 void gmx_sfree_impl(const char *name, const char *file, int line, T *ptr)
241 static_assert(std::is_pod<T>::value || std::is_void<T>::value,
242 "sfree() called on C++ type");
243 save_free(name, file, line, ptr);
245 /** C++ helper for sfree_aligned(). */
246 template <typename T> static inline
247 void gmx_sfree_aligned_impl(const char *name, const char *file, int line, T *ptr)
249 static_assert(std::is_pod<T>::value || std::is_void<T>::value,
250 "sfree_aligned() called on C++ type");
251 save_free_aligned(name, file, line, ptr);
253 /*! \} */
254 /*! \endcond */
256 /*! \def snew
257 * \brief
258 * Allocates memory for a given number of elements.
260 * \param[out] ptr Pointer to allocate.
261 * \param[in] nelem Number of elements to allocate.
263 * Allocates memory for \p nelem elements of type \p *ptr and sets this to
264 * \p ptr. The allocated memory is initialized to zeros.
266 * \hideinitializer
268 /*! \def srenew
269 * \brief
270 * Reallocates memory for a given number of elements.
272 * \param[in,out] ptr Pointer to allocate/reallocate.
273 * \param[in] nelem Number of elements to allocate.
275 * (Re)allocates memory for \p ptr such that it can hold \p nelem elements of
276 * type \p *ptr, and sets the new pointer to \p ptr.
277 * If \p ptr is `NULL`, memory is allocated as if it was new.
278 * If \p nelem is zero, \p ptr is freed (if not `NULL`).
279 * Note that the allocated memory is not initialized, unlike with snew().
281 * \hideinitializer
283 /*! \def smalloc
284 * \brief
285 * Allocates memory for a given number of bytes.
287 * \param[out] ptr Pointer to allocate.
288 * \param[in] size Number of bytes to allocate.
290 * Allocates memory for \p size bytes and sets this to \p ptr.
291 * The allocated memory is initialized to zero.
293 * \hideinitializer
295 /*! \def snew_aligned
296 * \brief
297 * Allocates aligned memory for a given number of elements.
299 * \param[out] ptr Pointer to allocate.
300 * \param[in] nelem Number of elements to allocate.
301 * \param[in] alignment Requested alignment in bytes.
303 * Allocates memory for \p nelem elements of type \p *ptr and sets this to
304 * \p ptr. The returned pointer is `alignment`-byte aligned.
305 * The allocated memory is initialized to zeros.
307 * The returned pointer should only be freed with sfree_aligned().
309 * \hideinitializer
311 /*! \def sfree
312 * \brief
313 * Frees memory referenced by \p ptr.
315 * \p ptr is allowed to be NULL, in which case nothing is done.
317 * \hideinitializer
319 /*! \def sfree_aligned
320 * \brief
321 * Frees aligned memory referenced by \p ptr.
323 * This must only be called with a pointer obtained through snew_aligned().
324 * \p ptr is allowed to be NULL, in which case nothing is done.
326 * \hideinitializer
329 /* C++ implementation */
330 #define snew(ptr, nelem) \
331 gmx_snew_impl(#ptr, __FILE__, __LINE__, (ptr), (nelem))
332 #define srenew(ptr, nelem) \
333 gmx_srenew_impl(#ptr, __FILE__, __LINE__, (ptr), (nelem))
334 #define smalloc(ptr, size) \
335 gmx_smalloc_impl(#ptr, __FILE__, __LINE__, (ptr), (size))
336 #define snew_aligned(ptr, nelem, alignment) \
337 gmx_snew_aligned_impl(#ptr, __FILE__, __LINE__, (ptr), (nelem), alignment)
338 #define sfree(ptr) \
339 gmx_sfree_impl(#ptr, __FILE__, __LINE__, (ptr))
340 #define sfree_aligned(ptr) \
341 gmx_sfree_aligned_impl(#ptr, __FILE__, __LINE__, (ptr))
343 /*! \brief
344 * Over allocation factor for memory allocations.
346 * Memory (re)allocation can be VERY slow, especially with some
347 * MPI libraries that replace the standard malloc and realloc calls.
348 * To avoid slow memory allocation we use over_alloc to set the memory
349 * allocation size for large data blocks. Since this scales the size
350 * with a factor, we use log(n) realloc calls instead of n.
351 * This can reduce allocation times from minutes to seconds.
353 * This factor leads to 4 realloc calls to double the array size.
355 constexpr float OVER_ALLOC_FAC = 1.19;
357 /*! \brief
358 * Turns over allocation for variable size atoms/cg/top arrays on or off,
359 * default is off.
361 * \todo
362 * This is mdrun-specific, so it might be better to put this and
363 * over_alloc_dd() much higher up.
365 void set_over_alloc_dd(gmx_bool set);
367 /*! \brief
368 * Returns new allocation count for domain decomposition allocations.
370 * Returns n when domain decomposition over allocation is off.
371 * Returns OVER_ALLOC_FAC*n + 100 when over allocation in on.
372 * This is to avoid frequent reallocation during domain decomposition in mdrun.
374 int over_alloc_dd(int n);
376 /** Over allocation for small data types: int, real etc. */
377 template<typename T>
378 constexpr T over_alloc_small(T n) { return OVER_ALLOC_FAC*n + 8000; }
380 /** Over allocation for large data types: complex structs */
381 template<typename T>
382 constexpr T over_alloc_large(T n) { return OVER_ALLOC_FAC*n + 1000; }
384 #endif