Gromacs  2026.0-dev-20241106-9ba7f4d
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Classes | Macros | Enumerations | Functions | Variables
atomdata.h File Reference
#include <cstdint>
#include <cstdio>
#include <memory>
#include <optional>
#include <vector>
#include "gromacs/gpu_utils/devicebuffer_datatype.h"
#include "gromacs/gpu_utils/hostallocator.h"
#include "gromacs/math/vectypes.h"
#include "gromacs/mdtypes/locality.h"
#include "gromacs/nbnxm/nbnxm_enums.h"
#include "gromacs/utility/arrayref.h"
#include "gromacs/utility/bitmask.h"
#include "gromacs/utility/real.h"
+ Include dependency graph for atomdata.h:
+ This graph shows which files directly or indirectly include this file:

Description

Functionality for per-atom data in the nbnxm module.

Author
Berk Hess hess@.nosp@m.kth..nosp@m.se

Classes

class  gmx::EnergyAccumulator< bool, bool >
 Base energy accumulator class, only specializations are used. More...
 
struct  gmx::nbnxn_atomdata_output_t
 Struct that holds force and energy output buffers. More...
 
struct  gmx::nbnxn_atomdata_t
 Struct that stores atom related data for the nbnxn module. More...
 
struct  gmx::nbnxn_atomdata_t::Params
 The actual atom data parameter values. More...
 
struct  gmx::nbnxn_atomdata_t::SimdMasks
 Diagonal and topology exclusion helper data for all SIMD kernels. More...
 

Macros

#define NBNXN_BUFFERFLAG_SIZE   16
 Block size in atoms for the non-bonded thread force-buffer reduction. More...
 
#define NBNXN_BUFFERFLAG_MAX_THREADS   (BITMASK_SIZE)
 We store the reduction flags as gmx_bitmask_t. This limits the number of flags to BITMASK_SIZE.
 

Enumerations

enum  { nbatXYZ, nbatXYZQ, nbatX4, nbatX8 }
 
enum  gmx::LJCombinationRule : int { gmx::LJCombinationRule::Geometric, gmx::LJCombinationRule::LorentzBerthelot, gmx::LJCombinationRule::None, gmx::LJCombinationRule::Count }
 LJ combination rules. More...
 

Functions

template<int packSize>
static int gmx::atom_to_x_index (int a)
 Returns the index in a coordinate array corresponding to atom a.
 
const char * gmx::enumValueToString (LJCombinationRule enumValue)
 String corresponding to LJ combination rule.
 
void gmx::copy_rvec_to_nbat_real (const int *a, int na, int na_round, const rvec *x, int nbatFormat, real *xnb, int a0)
 Copy na rvec elements from x to xnb using nbatFormat, start dest a0, and fills up to na_round with coordinates that are far away.
 
void gmx::nbnxn_atomdata_set (nbnxn_atomdata_t *nbat, const GridSet &gridSet, ArrayRef< const int > atomTypes, ArrayRef< const real > atomCharges, ArrayRef< const int32_t > atomInfo)
 Sets the atomdata after pair search.
 
void gmx::nbnxn_atomdata_copy_shiftvec (bool dynamic_box, ArrayRef< RVec > shift_vec, nbnxn_atomdata_t *nbat)
 Copy the shift vectors to nbat.
 
void gmx::nbnxn_atomdata_copy_x_to_nbat_x (const GridSet &gridSet, AtomLocality locality, const rvec *coordinates, nbnxn_atomdata_t *nbat)
 Transform coordinates to xbat layout. More...
 
void gmx::nbnxn_atomdata_x_to_nbat_x_gpu (const GridSet &gridSet, AtomLocality locality, NbnxmGpu *gpu_nbv, DeviceBuffer< RVec > d_x, GpuEventSynchronizer *xReadyOnDevice)
 Transform coordinates to xbat layout on GPU. More...
 
void gmx::nbnxn_atomdata_add_nbat_fshift_to_fshift (const nbnxn_atomdata_t &nbat, ArrayRef< RVec > fshift)
 Add the fshift force stored in nbat to fshift.
 

Variables

static constexpr int gmx::STRIDE_XYZ = 3
 Stride for coordinate/force arrays with xyz coordinate storage.
 
static constexpr int gmx::STRIDE_XYZQ = 4
 Stride for coordinate/force arrays with xyzq coordinate storage.
 
static constexpr int gmx::c_packX4 = 4
 Size of packs of x, y or z with SIMD 4-grouped packed coordinates/forces.
 
static constexpr int gmx::c_packX8 = 8
 Size of packs of x, y or z with SIMD 8-grouped packed coordinates/forces.
 
static constexpr int gmx::STRIDE_P4 = DIM * c_packX4
 Stridefor a pack of 4 coordinates/forces.
 
static constexpr int gmx::STRIDE_P8 = DIM * c_packX8
 Stridefor a pack of 8 coordinates/forces.
 

Macro Definition Documentation

#define NBNXN_BUFFERFLAG_SIZE   16

Block size in atoms for the non-bonded thread force-buffer reduction.

Should be a multiple of all cell and x86 SIMD sizes (i.e. 2, 4 and 8). Should be small to reduce the reduction and zeroing cost, but too small will result in overhead. Currently the block size is NBNXN_BUFFERFLAG_SIZE*3*sizeof(real)=192 bytes.