Gromacs
2020.4
|
C-style memory allocation routines for GROMACS.
This header provides macros snew(), srenew(), smalloc(), and sfree() for C-style memory management. Additionally, snew_aligned() and sfree_aligned() are provided for managing memory with a specified byte alignment.
If an allocation fails, the program is halted by calling gmx_fatal(), which outputs source file and line number and the name of the variable involved. This frees calling code from the trouble of checking the result of the allocations everywhere. It also provides a location for centrally logging memory allocations for diagnosing memory usage (currently can only enabled by changing the source code). Additionally, sfree() works also with a NULL
parameter, which standard free() does not.
The macros forward the calls to functions save_malloc(), save_calloc(), save_realloc(), save_free(), save_calloc_aligned(), and save_free_aligned(). There are a few low-level locations in GROMACS that call these directly, but generally the macros should be used. save_malloc_aligned() exists for this purpose, although there is no macro to invoke it.
Macros | |
#define | snew(ptr, nelem) |
Allocates memory for a given number of elements. More... | |
#define | srenew(ptr, nelem) |
Reallocates memory for a given number of elements. More... | |
#define | smalloc(ptr, size) |
Allocates memory for a given number of bytes. More... | |
#define | snew_aligned(ptr, nelem, alignment) |
Allocates aligned memory for a given number of elements. More... | |
#define | sfree(ptr) |
Frees memory referenced by ptr . More... | |
#define | sfree_aligned(ptr) |
Frees aligned memory referenced by ptr . More... | |
Functions | |
void * | save_malloc (const char *name, const char *file, int line, size_t size) |
GROMACS wrapper for malloc(). More... | |
void * | save_calloc (const char *name, const char *file, int line, size_t nelem, size_t elsize) |
GROMACS wrapper for calloc(). More... | |
void * | save_realloc (const char *name, const char *file, int line, void *ptr, size_t nelem, size_t elsize) |
GROMACS wrapper for realloc(). More... | |
void | save_free (const char *name, const char *file, int line, void *ptr) |
GROMACS wrapper for free(). More... | |
void * | save_malloc_aligned (const char *name, const char *file, int line, size_t nelem, size_t elsize, size_t alignment) |
GROMACS wrapper for allocating aligned memory. More... | |
void * | save_calloc_aligned (const char *name, const char *file, int line, size_t nelem, size_t elsize, size_t alignment) |
GROMACS wrapper for allocating zero-initialized aligned memory. More... | |
void | save_free_aligned (const char *name, const char *file, int line, void *ptr) |
GROMACS wrapper for freeing aligned memory. More... | |
void | set_over_alloc_dd (gmx_bool set) |
Turns over allocation for variable size atoms/cg/top arrays on or off, default is off. More... | |
int | over_alloc_dd (int n) |
Returns new allocation count for domain decomposition allocations. More... | |
template<typename T > | |
constexpr T | over_alloc_small (T n) |
Over allocation for small data types: int, real etc. More... | |
template<typename T > | |
constexpr T | over_alloc_large (T n) |
Over allocation for large data types: complex structs. | |
Implementation templates for C++ memory allocation macros | |
These templates are used to implement the snew() etc. macros for C++, where an explicit cast is needed from Having these as The size cannot be passed as a parameter to the function either, since that provokes warnings from cppcheck for some invocations, where a complex expression is passed as | |
template<typename T > | |
static void | gmx_snew_impl (const char *name, const char *file, int line, T *&ptr, size_t nelem) |
C++ helper for snew(). More... | |
template<typename T > | |
static void | gmx_srenew_impl (const char *name, const char *file, int line, T *&ptr, size_t nelem) |
C++ helper for srenew(). More... | |
template<typename T > | |
static void | gmx_smalloc_impl (const char *name, const char *file, int line, T *&ptr, size_t size) |
C++ helper for smalloc(). More... | |
template<typename T > | |
static void | gmx_snew_aligned_impl (const char *name, const char *file, int line, T *&ptr, size_t nelem, size_t alignment) |
C++ helper for snew_aligned(). More... | |
template<typename T > | |
static void | gmx_sfree_impl (const char *name, const char *file, int line, T *ptr) |
C++ helper for sfree(). More... | |
template<typename T > | |
static void | gmx_sfree_aligned_impl (const char *name, const char *file, int line, T *ptr) |
C++ helper for sfree_aligned(). More... | |
Variables | |
constexpr float | OVER_ALLOC_FAC = 1.19F |
Over allocation factor for memory allocations. More... | |
#define sfree | ( | ptr | ) |
Frees memory referenced by ptr
.
ptr
is allowed to be NULL, in which case nothing is done.
#define sfree_aligned | ( | ptr | ) |
Frees aligned memory referenced by ptr
.
This must only be called with a pointer obtained through snew_aligned(). ptr
is allowed to be NULL, in which case nothing is done.
#define smalloc | ( | ptr, | |
size | |||
) |
Allocates memory for a given number of bytes.
[out] | ptr | Pointer to allocate. |
[in] | size | Number of bytes to allocate. |
Allocates memory for size
bytes and sets this to ptr
. The allocated memory is initialized to zero.
#define snew | ( | ptr, | |
nelem | |||
) |
Allocates memory for a given number of elements.
[out] | ptr | Pointer to allocate. |
[in] | nelem | Number of elements to allocate. |
Allocates memory for nelem
elements of type *ptr
and sets this to ptr
. The allocated memory is initialized to zeros.
#define snew_aligned | ( | ptr, | |
nelem, | |||
alignment | |||
) |
Allocates aligned memory for a given number of elements.
[out] | ptr | Pointer to allocate. |
[in] | nelem | Number of elements to allocate. |
[in] | alignment | Requested alignment in bytes. |
Allocates memory for nelem
elements of type *ptr
and sets this to ptr
. The returned pointer is alignment
-byte aligned. The allocated memory is initialized to zeros.
The returned pointer should only be freed with sfree_aligned().
#define srenew | ( | ptr, | |
nelem | |||
) |
Reallocates memory for a given number of elements.
[in,out] | ptr | Pointer to allocate/reallocate. |
[in] | nelem | Number of elements to allocate. |
(Re)allocates memory for ptr
such that it can hold nelem
elements of type *ptr
, and sets the new pointer to ptr
. If ptr
is NULL
, memory is allocated as if it was new. If nelem
is zero, ptr
is freed (if not NULL
). Note that the allocated memory is not initialized, unlike with snew().
|
inlinestatic |
C++ helper for sfree_aligned().
|
inlinestatic |
C++ helper for sfree().
|
inlinestatic |
C++ helper for smalloc().
|
inlinestatic |
C++ helper for snew_aligned().
|
inlinestatic |
C++ helper for snew().
|
inlinestatic |
C++ helper for srenew().
int over_alloc_dd | ( | int | n | ) |
Returns new allocation count for domain decomposition allocations.
Returns n when domain decomposition over allocation is off. Returns OVER_ALLOC_FAC*n + 100 when over allocation in on. This is to avoid frequent reallocation during domain decomposition in mdrun.
constexpr T over_alloc_small | ( | T | n | ) |
Over allocation for small data types: int, real etc.
void* save_calloc | ( | const char * | name, |
const char * | file, | ||
int | line, | ||
size_t | nelem, | ||
size_t | elsize | ||
) |
GROMACS wrapper for calloc().
[in] | name | Variable name identifying the allocation. |
[in] | file | Source code file where the allocation originates from. |
[in] | line | Source code line where the allocation originates from. |
[in] | nelem | Number of elements to allocate. |
[in] | elsize | Number of bytes per element. |
This should generally be called through snew(), not directly.
void* save_calloc_aligned | ( | const char * | name, |
const char * | file, | ||
int | line, | ||
size_t | nelem, | ||
size_t | elsize, | ||
size_t | alignment | ||
) |
GROMACS wrapper for allocating zero-initialized aligned memory.
[in] | name | Variable name identifying the allocation. |
[in] | file | Source code file where the allocation originates from. |
[in] | line | Source code line where the allocation originates from. |
[in] | nelem | Number of elements to allocate. |
[in] | elsize | Number of bytes per element. |
[in] | alignment | Requested alignment in bytes. |
alignment
-byte boundary.This should generally be called through snew_aligned(), not directly.
The returned pointer should only be freed with a call to save_free_aligned().
void save_free | ( | const char * | name, |
const char * | file, | ||
int | line, | ||
void * | ptr | ||
) |
GROMACS wrapper for free().
[in] | name | Variable name identifying the deallocation. |
[in] | file | Source code file where the deallocation originates from. |
[in] | line | Source code line where the deallocation originates from. |
[in] | ptr | Pointer to the allocated memory (can be NULL). |
If ptr
is NULL, does nothing. This should generally be called through sfree(), not directly. This never fails.
void save_free_aligned | ( | const char * | name, |
const char * | file, | ||
int | line, | ||
void * | ptr | ||
) |
GROMACS wrapper for freeing aligned memory.
[in] | name | Variable name identifying the deallocation. |
[in] | file | Source code file where the deallocation originates from. |
[in] | line | Source code line where the deallocation originates from. |
[in] | ptr | Pointer to the allocated memory (can be NULL). |
If ptr
is NULL, does nothing. ptr
should have been allocated with save_malloc_aligned() or save_calloc_aligned(). This should generally be called through sfree_aligned(), not directly. This never fails.
void* save_malloc | ( | const char * | name, |
const char * | file, | ||
int | line, | ||
size_t | size | ||
) |
GROMACS wrapper for malloc().
[in] | name | Variable name identifying the allocation. |
[in] | file | Source code file where the allocation originates from. |
[in] | line | Source code line where the allocation originates from. |
[in] | size | Number of bytes to allocate. |
This should generally be called through smalloc(), not directly.
void* save_malloc_aligned | ( | const char * | name, |
const char * | file, | ||
int | line, | ||
size_t | nelem, | ||
size_t | elsize, | ||
size_t | alignment | ||
) |
GROMACS wrapper for allocating aligned memory.
[in] | name | Variable name identifying the allocation. |
[in] | file | Source code file where the allocation originates from. |
[in] | line | Source code line where the allocation originates from. |
[in] | nelem | Number of elements to allocate. |
[in] | elsize | Number of bytes per element. |
[in] | alignment | Requested alignment in bytes. |
alignment
-byte boundary.There is no macro that invokes this function.
The returned pointer should only be freed with a call to save_free_aligned().
void* save_realloc | ( | const char * | name, |
const char * | file, | ||
int | line, | ||
void * | ptr, | ||
size_t | nelem, | ||
size_t | elsize | ||
) |
GROMACS wrapper for realloc().
[in] | name | Variable name identifying the allocation. |
[in] | file | Source code file where the allocation originates from. |
[in] | line | Source code line where the allocation originates from. |
[in] | ptr | Pointer to the previously allocated memory (can be NULL). |
[in] | nelem | Number of elements to allocate. |
[in] | elsize | Number of bytes per element. |
As with realloc(), if ptr
is NULL, memory is allocated as if malloc() was called. This should generally be called through srenew(), not directly.
Note that the allocated memory is not initialized to zero.
void set_over_alloc_dd | ( | gmx_bool | set | ) |
Turns over allocation for variable size atoms/cg/top arrays on or off, default is off.
constexpr float OVER_ALLOC_FAC = 1.19F |
Over allocation factor for memory allocations.
Memory (re)allocation can be VERY slow, especially with some MPI libraries that replace the standard malloc and realloc calls. To avoid slow memory allocation we use over_alloc to set the memory allocation size for large data blocks. Since this scales the size with a factor, we use log(n) realloc calls instead of n. This can reduce allocation times from minutes to seconds.
This factor leads to 4 realloc calls to double the array size.