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