Gromacs  2024.4
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Macros | Functions | Variables
smalloc.h File Reference
#include <cstddef>
#include <type_traits>
+ Include dependency graph for smalloc.h:
+ This graph shows which files directly or indirectly include this file:

Description

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 (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.
 

Variables

constexpr float OVER_ALLOC_FAC = 1.19F
 Over allocation factor for memory allocations. More...
 

Macro Definition Documentation

#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.

Parameters
[out]ptrPointer to allocate.
[in]sizeNumber 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.

Parameters
[out]ptrPointer to allocate.
[in]nelemNumber 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.

Parameters
[out]ptrPointer to allocate.
[in]nelemNumber of elements to allocate.
[in]alignmentRequested 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.

Parameters
[in,out]ptrPointer to allocate/reallocate.
[in]nelemNumber 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().

Function Documentation

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.

template<typename T >
constexpr T over_alloc_small ( 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().

Parameters
[in]nameVariable name identifying the allocation.
[in]fileSource code file where the allocation originates from.
[in]lineSource code line where the allocation originates from.
[in]nelemNumber of elements to allocate.
[in]elsizeNumber of bytes per element.
Returns
Pointer to the allocated space.

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.

Parameters
[in]nameVariable name identifying the allocation.
[in]fileSource code file where the allocation originates from.
[in]lineSource code line where the allocation originates from.
[in]nelemNumber of elements to allocate.
[in]elsizeNumber of bytes per element.
[in]alignmentRequested alignment in bytes.
Returns
Pointer to the allocated space, aligned at 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().

Parameters
[in]nameVariable name identifying the deallocation.
[in]fileSource code file where the deallocation originates from.
[in]lineSource code line where the deallocation originates from.
[in]ptrPointer 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.

Parameters
[in]nameVariable name identifying the deallocation.
[in]fileSource code file where the deallocation originates from.
[in]lineSource code line where the deallocation originates from.
[in]ptrPointer 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().

Parameters
[in]nameVariable name identifying the allocation.
[in]fileSource code file where the allocation originates from.
[in]lineSource code line where the allocation originates from.
[in]sizeNumber of bytes to allocate.
Returns
Pointer to the allocated space.

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.

Parameters
[in]nameVariable name identifying the allocation.
[in]fileSource code file where the allocation originates from.
[in]lineSource code line where the allocation originates from.
[in]nelemNumber of elements to allocate.
[in]elsizeNumber of bytes per element.
[in]alignmentRequested alignment in bytes.
Returns
Pointer to the allocated space, aligned at 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().

Parameters
[in]nameVariable name identifying the allocation.
[in]fileSource code file where the allocation originates from.
[in]lineSource code line where the allocation originates from.
[in]ptrPointer to the previously allocated memory (can be NULL).
[in]nelemNumber of elements to allocate.
[in]elsizeNumber of bytes per element.
Returns
Pointer to the allocated space.

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 ( bool  set)

Turns over allocation for variable size atoms/cg/top arrays on or off, default is off.

Todo:
This is mdrun-specific, so it might be better to put this and over_alloc_dd() much higher up.

Variable Documentation

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.