1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 /*! \file
38  * \brief
39  * C-style memory allocation routines for \Gromacs.
40  *
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.
44  *
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.
52  *
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.
59  *
60  * \inpublicapi
61  * \ingroup module_utility
62  */
63 #ifndef GMX_UTILITY_SMALLOC_H
64 #define GMX_UTILITY_SMALLOC_H
65 
66 #include <stddef.h>
67 
68 #include "gromacs/utility/basedefinitions.h"
69 
70 /*! \brief
71  * \Gromacs wrapper for malloc().
72  *
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.
78  *
79  * This should generally be called through smalloc(), not directly.
80  */
81 void* save_malloc(const char* name, const char* file, int line, size_t size);
82 /*! \brief
83  * \Gromacs wrapper for calloc().
84  *
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.
91  *
92  * This should generally be called through snew(), not directly.
93  */
94 void* save_calloc(const char* name, const char* file, int line, size_t nelem, size_t elsize);
95 /*! \brief
96  * \Gromacs wrapper for realloc().
97  *
98  * \param[in] name   Variable name identifying the allocation.
99  * \param[in] file   Source code file where the allocation originates from.
100  * \param[in] line   Source code line where the allocation originates from.
101  * \param[in] ptr    Pointer to the previously allocated memory (can be NULL).
102  * \param[in] nelem  Number of elements to allocate.
103  * \param[in] elsize Number of bytes per element.
104  * \returns   Pointer to the allocated space.
105  *
106  * As with realloc(), if \p ptr is NULL, memory is allocated as if malloc() was
107  * called.
108  * This should generally be called through srenew(), not directly.
109  *
110  * Note that the allocated memory is not initialized to zero.
111  */
112 void* save_realloc(const char* name, const char* file, int line, void* ptr, size_t nelem, size_t elsize);
113 /*! \brief
114  * \Gromacs wrapper for free().
115  *
116  * \param[in] name   Variable name identifying the deallocation.
117  * \param[in] file   Source code file where the deallocation originates from.
118  * \param[in] line   Source code line where the deallocation originates from.
119  * \param[in] ptr    Pointer to the allocated memory (can be NULL).
120  *
121  * If \p ptr is NULL, does nothing.
122  * This should generally be called through sfree(), not directly.
123  * This never fails.
124  */
125 void save_free(const char* name, const char* file, int line, void* ptr);
126 
127 /*! \brief
128  * \Gromacs wrapper for allocating aligned memory.
129  *
130  * \param[in] name   Variable name identifying the allocation.
131  * \param[in] file   Source code file where the allocation originates from.
132  * \param[in] line   Source code line where the allocation originates from.
133  * \param[in] nelem  Number of elements to allocate.
134  * \param[in] elsize Number of bytes per element.
135  * \param[in] alignment Requested alignment in bytes.
136  * \returns   Pointer to the allocated space, aligned at `alignment`-byte
137  *     boundary.
138  *
139  * There is no macro that invokes this function.
140  *
141  * The returned pointer should only be freed with a call to save_free_aligned().
142  */
143 void* save_malloc_aligned(const char* name, const char* file, int line, size_t nelem, size_t elsize, size_t alignment);
144 /*! \brief
145  * \Gromacs wrapper for allocating zero-initialized aligned memory.
146  *
147  * \param[in] name   Variable name identifying the allocation.
148  * \param[in] file   Source code file where the allocation originates from.
149  * \param[in] line   Source code line where the allocation originates from.
150  * \param[in] nelem  Number of elements to allocate.
151  * \param[in] elsize Number of bytes per element.
152  * \param[in] alignment Requested alignment in bytes.
153  * \returns   Pointer to the allocated space, aligned at `alignment`-byte
154  *     boundary.
155  *
156  * This should generally be called through snew_aligned(), not directly.
157  *
158  * The returned pointer should only be freed with a call to save_free_aligned().
159  */
160 void* save_calloc_aligned(const char* name, const char* file, int line, size_t nelem, size_t elsize, size_t alignment);
161 /*! \brief
162  * \Gromacs wrapper for freeing aligned memory.
163  *
164  * \param[in] name   Variable name identifying the deallocation.
165  * \param[in] file   Source code file where the deallocation originates from.
166  * \param[in] line   Source code line where the deallocation originates from.
167  * \param[in] ptr    Pointer to the allocated memory (can be NULL).
168  *
169  * If \p ptr is NULL, does nothing.
170  * \p ptr should have been allocated with save_malloc_aligned() or
171  * save_calloc_aligned().
172  * This should generally be called through sfree_aligned(), not directly.
173  * This never fails.
174  */
175 void save_free_aligned(const char* name, const char* file, int line, void* ptr);
176 
177 #include <type_traits>
178 
179 /*! \cond internal */
180 /*! \name Implementation templates for C++ memory allocation macros
181  *
182  * These templates are used to implement the snew() etc. macros for C++, where
183  * an explicit cast is needed from `void *` (the return value of the allocation
184  * wrapper functions) to the type of \p ptr.
185  *
186  * Having these as `static` avoid some obscure bugs if several files define
187  * distinct data structures with identical names and allocate memory for them
188  * using snew().  By the C++ standard, such declarations cause undefined
189  * behavior, but can be difficult to spot in the existing C code.
190  * Without the `static` (and if the compiler does not inline the calls), the
191  * linker cannot that data structures with identical names are actually
192  * different and links calls to these template functions incorrectly, which can
193  * result in allocation of an incorrect amount of memory if the element size is
194  * computed within the function.
195  *
196  * The size cannot be passed as a parameter to the function either, since that
197  * provokes warnings from cppcheck for some invocations, where a complex
198  * expression is passed as \p ptr.
199  */
200 /*! \{ */
201 /** C++ helper for snew(). */
202 template<typename T>
gmx_snew_impl(const char * name,const char * file,int line,T * & ptr,size_t nelem)203 static inline void gmx_snew_impl(const char* name, const char* file, int line, T*& ptr, size_t nelem)
204 {
205     // TODO: Use std::is_pod_v when CUDA 11 is a requirement.
206     static_assert(std::is_pod<T>::value, "snew() called on C++ type");
207     // NOLINTNEXTLINE bugprone-sizeof-expression
208     ptr = static_cast<T*>(save_calloc(name, file, line, nelem, sizeof(T)));
209 }
210 /** C++ helper for srenew(). */
211 template<typename T>
gmx_srenew_impl(const char * name,const char * file,int line,T * & ptr,size_t nelem)212 static inline void gmx_srenew_impl(const char* name, const char* file, int line, T*& ptr, size_t nelem)
213 {
214     // TODO: Use std::is_pod_v when CUDA 11 is a requirement.
215     static_assert(std::is_pod<T>::value, "srenew() called on C++ type");
216     // NOLINTNEXTLINE bugprone-sizeof-expression
217     ptr = static_cast<T*>(save_realloc(name, file, line, ptr, nelem, sizeof(T)));
218 }
219 /** C++ helper for smalloc(). */
220 template<typename T>
gmx_smalloc_impl(const char * name,const char * file,int line,T * & ptr,size_t size)221 static inline void gmx_smalloc_impl(const char* name, const char* file, int line, T*& ptr, size_t size)
222 {
223     // TODO: Use std::is_pod_v when CUDA 11 is a requirement.
224     static_assert(std::is_pod<T>::value, "smalloc() called on C++ type");
225     ptr = static_cast<T*>(save_malloc(name, file, line, size));
226 }
227 /** C++ helper for snew_aligned(). */
228 template<typename T>
229 static inline void
gmx_snew_aligned_impl(const char * name,const char * file,int line,T * & ptr,size_t nelem,size_t alignment)230 gmx_snew_aligned_impl(const char* name, const char* file, int line, T*& ptr, size_t nelem, size_t alignment)
231 {
232     // TODO: Use std::is_pod_v when CUDA 11 is a requirement.
233     static_assert(std::is_pod<T>::value, "snew_aligned() called on C++ type");
234     ptr = static_cast<T*>(save_calloc_aligned(name, file, line, nelem, sizeof(T), alignment));
235 }
236 /** C++ helper for sfree(). */
237 template<typename T>
gmx_sfree_impl(const char * name,const char * file,int line,T * ptr)238 static inline void gmx_sfree_impl(const char* name, const char* file, int line, T* ptr)
239 {
240     // TODO: Use std::is_pod_v and std::is_void_v when CUDA 11 is a requirement.
241     static_assert(std::is_pod<T>::value || std::is_void<T>::value, "sfree() called on C++ type");
242     save_free(name, file, line, ptr);
243 }
244 /** C++ helper for sfree_aligned(). */
245 template<typename T>
gmx_sfree_aligned_impl(const char * name,const char * file,int line,T * ptr)246 static inline void gmx_sfree_aligned_impl(const char* name, const char* file, int line, T* ptr)
247 {
248     // TODO: Use std::is_pod_v and std::is_void_v when CUDA 11 is a requirement.
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);
252 }
253 /*! \} */
254 /*! \endcond */
255 
256 /*! \def snew
257  * \brief
258  * Allocates memory for a given number of elements.
259  *
260  * \param[out] ptr   Pointer to allocate.
261  * \param[in]  nelem Number of elements to allocate.
262  *
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.
265  *
266  * \hideinitializer
267  */
268 /*! \def srenew
269  * \brief
270  * Reallocates memory for a given number of elements.
271  *
272  * \param[in,out] ptr   Pointer to allocate/reallocate.
273  * \param[in]     nelem Number of elements to allocate.
274  *
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().
280  *
281  * \hideinitializer
282  */
283 /*! \def smalloc
284  * \brief
285  * Allocates memory for a given number of bytes.
286  *
287  * \param[out] ptr  Pointer to allocate.
288  * \param[in]  size Number of bytes to allocate.
289  *
290  * Allocates memory for \p size bytes and sets this to \p ptr.
291  * The allocated memory is initialized to zero.
292  *
293  * \hideinitializer
294  */
295 /*! \def snew_aligned
296  * \brief
297  * Allocates aligned memory for a given number of elements.
298  *
299  * \param[out] ptr       Pointer to allocate.
300  * \param[in]  nelem     Number of elements to allocate.
301  * \param[in]  alignment Requested alignment in bytes.
302  *
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.
306  *
307  * The returned pointer should only be freed with sfree_aligned().
308  *
309  * \hideinitializer
310  */
311 /*! \def sfree
312  * \brief
313  * Frees memory referenced by \p ptr.
314  *
315  * \p ptr is allowed to be NULL, in which case nothing is done.
316  *
317  * \hideinitializer
318  */
319 /*! \def sfree_aligned
320  * \brief
321  * Frees aligned memory referenced by \p ptr.
322  *
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.
325  *
326  * \hideinitializer
327  */
328 
329 /* C++ implementation */
330 #define snew(ptr, nelem) gmx_snew_impl(#ptr, __FILE__, __LINE__, (ptr), (nelem))
331 #define srenew(ptr, nelem) gmx_srenew_impl(#ptr, __FILE__, __LINE__, (ptr), (nelem))
332 #define smalloc(ptr, size) gmx_smalloc_impl(#ptr, __FILE__, __LINE__, (ptr), (size))
333 #define snew_aligned(ptr, nelem, alignment) \
334     gmx_snew_aligned_impl(#ptr, __FILE__, __LINE__, (ptr), (nelem), alignment)
335 #define sfree(ptr) gmx_sfree_impl(#ptr, __FILE__, __LINE__, (ptr))
336 #define sfree_aligned(ptr) gmx_sfree_aligned_impl(#ptr, __FILE__, __LINE__, (ptr))
337 
338 /*! \brief
339  * Over allocation factor for memory allocations.
340  *
341  * Memory (re)allocation can be VERY slow, especially with some
342  * MPI libraries that replace the standard malloc and realloc calls.
343  * To avoid slow memory allocation we use over_alloc to set the memory
344  * allocation size for large data blocks. Since this scales the size
345  * with a factor, we use log(n) realloc calls instead of n.
346  * This can reduce allocation times from minutes to seconds.
347  *
348  * This factor leads to 4 realloc calls to double the array size.
349  */
350 constexpr float OVER_ALLOC_FAC = 1.19F;
351 
352 /*! \brief
353  * Turns over allocation for variable size atoms/cg/top arrays on or off,
354  * default is off.
355  *
356  * \todo
357  * This is mdrun-specific, so it might be better to put this and
358  * over_alloc_dd() much higher up.
359  */
360 void set_over_alloc_dd(gmx_bool set);
361 
362 /*! \brief
363  * Returns new allocation count for domain decomposition allocations.
364  *
365  * Returns n when domain decomposition over allocation is off.
366  * Returns OVER_ALLOC_FAC*n + 100 when over allocation in on.
367  * This is to avoid frequent reallocation during domain decomposition in mdrun.
368  */
369 int over_alloc_dd(int n);
370 
371 /** Over allocation for small data types: int, real etc. */
372 template<typename T>
over_alloc_small(T n)373 constexpr T over_alloc_small(T n)
374 {
375     return OVER_ALLOC_FAC * n + 8000;
376 }
377 
378 /** Over allocation for large data types: complex structs */
379 template<typename T>
over_alloc_large(T n)380 constexpr T over_alloc_large(T n)
381 {
382     return OVER_ALLOC_FAC * n + 1000;
383 }
384 
385 #endif
386