Gromacs  2023.2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members | Public Member Functions | Public Attributes
nonbonded_verlet_t Struct Reference

#include <gromacs/nbnxm/nbnxm.h>

Description

Top-level non-bonded data structure for the Verlet-type cut-off scheme.

Public Member Functions

 nonbonded_verlet_t (std::unique_ptr< PairlistSets > pairlistSets, std::unique_ptr< PairSearch > pairSearch, std::unique_ptr< nbnxn_atomdata_t > nbat, const Nbnxm::KernelSetup &kernelSetup, std::unique_ptr< ExclusionChecker > exclusionChecker, NbnxmGpu *gpu_nbv, gmx_wallcycle *wcycle)
 Constructs an object from its components. More...
 
 nonbonded_verlet_t (std::unique_ptr< PairlistSets > pairlistSets, std::unique_ptr< PairSearch > pairSearch, std::unique_ptr< nbnxn_atomdata_t > nbat, const Nbnxm::KernelSetup &kernelSetup, NbnxmGpu *gpu_nbv)
 Constructs an object from its, minimal, components. More...
 
bool useGpu () const
 Returns whether a GPU is use for the non-bonded calculations.
 
bool emulateGpu () const
 Returns whether a GPU is emulated for the non-bonded calculations.
 
bool pairlistIsSimple () const
 Return whether the pairlist is of simple, CPU type.
 
gmx::ArrayRef< const int > getLocalAtomOrder () const
 Returns the order of the local atoms on the grid.
 
void setLocalAtomOrder () const
 Sets the order of the local atoms to the order grid atom ordering.
 
gmx::ArrayRef< const int > getGridIndices () const
 Returns the index position of the atoms on the search grid.
 
void constructPairlist (gmx::InteractionLocality iLocality, const gmx::ListOfLists< int > &exclusions, int64_t step, t_nrnb *nrnb) const
 Constructs the pairlist for the given locality. More...
 
void setAtomProperties (gmx::ArrayRef< const int > atomTypes, gmx::ArrayRef< const real > atomCharges, gmx::ArrayRef< const int64_t > atomInfo) const
 Updates all the atom properties in Nbnxm.
 
void convertCoordinates (gmx::AtomLocality locality, gmx::ArrayRef< const gmx::RVec > coordinates)
 Convert the coordinates to NBNXM format for the given locality. More...
 
void convertCoordinatesGpu (gmx::AtomLocality locality, DeviceBuffer< gmx::RVec > d_x, GpuEventSynchronizer *xReadyOnDevice)
 Convert the coordinates to NBNXM format on the GPU for the given locality. More...
 
void atomdata_init_copy_x_to_nbat_x_gpu () const
 Init for GPU version of setup coordinates in Nbnxm.
 
const PairlistSetspairlistSets () const
 Returns a reference to the pairlist sets.
 
bool isDynamicPruningStepCpu (int64_t step) const
 Returns whether step is a dynamic list pruning step, for CPU lists.
 
bool isDynamicPruningStepGpu (int64_t step) const
 Returns whether step is a dynamic list pruning step, for GPU lists.
 
void dispatchPruneKernelCpu (gmx::InteractionLocality iLocality, gmx::ArrayRef< const gmx::RVec > shift_vec) const
 Dispatches the dynamic pruning kernel for the given locality, for CPU lists.
 
void dispatchPruneKernelGpu (int64_t step)
 Dispatches the dynamic pruning kernel for GPU lists.
 
void dispatchNonbondedKernel (gmx::InteractionLocality iLocality, const interaction_const_t &ic, const gmx::StepWorkload &stepWork, int clearF, gmx::ArrayRef< const gmx::RVec > shiftvec, gmx::ArrayRef< real > repulsionDispersionSR, gmx::ArrayRef< real > CoulombSR, t_nrnb *nrnb) const
 Executes the non-bonded kernel of the GPU or launches it on the GPU.
 
void dispatchFreeEnergyKernels (const gmx::ArrayRefWithPadding< const gmx::RVec > &coords, gmx::ForceWithShiftForces *forceWithShiftForces, bool useSimd, int ntype, const interaction_const_t &ic, gmx::ArrayRef< const gmx::RVec > shiftvec, gmx::ArrayRef< const real > nbfp, gmx::ArrayRef< const real > nbfp_grid, gmx::ArrayRef< const real > chargeA, gmx::ArrayRef< const real > chargeB, gmx::ArrayRef< const int > typeA, gmx::ArrayRef< const int > typeB, t_lambda *fepvals, gmx::ArrayRef< const real > lambda, gmx_enerdata_t *enerd, const gmx::StepWorkload &stepWork, t_nrnb *nrnb)
 Executes the non-bonded free-energy kernels, local + non-local, always runs on the CPU.
 
void atomdata_add_nbat_f_to_f (gmx::AtomLocality locality, gmx::ArrayRef< gmx::RVec > force)
 Add the forces stored in nbat to f, zeros the forces in nbat. More...
 
int getNumAtoms (gmx::AtomLocality locality) const
 Get the number of atoms for a given locality. More...
 
const Nbnxm::KernelSetup & kernelSetup () const
 Return the kernel setup.
 
real pairlistInnerRadius () const
 Returns the outer radius for the pair list.
 
real pairlistOuterRadius () const
 Returns the outer radius for the pair list.
 
void changePairlistRadii (real rlistOuter, real rlistInner) const
 Changes the pair-list outer and inner radius.
 
void setupGpuShortRangeWork (const gmx::ListedForcesGpu *listedForcesGpu, gmx::InteractionLocality iLocality) const
 Set up internal flags that indicate what type of short-range work there is.
 
void setupFepThreadedForceBuffer (int numAtomsForce)
 

Public Attributes

std::unique_ptr< PairlistSetspairlistSets_
 All data related to the pair lists.
 
std::unique_ptr< PairSearchpairSearch_
 Working data for constructing the pairlists.
 
std::unique_ptr< nbnxn_atomdata_t > nbat
 Atom data.
 
NbnxmGpu * gpu_nbv
 GPU Nbnxm data, only used with a physical GPU (TODO: use unique_ptr)
 

Constructor & Destructor Documentation

nonbonded_verlet_t::nonbonded_verlet_t ( std::unique_ptr< PairlistSets pairlistSets,
std::unique_ptr< PairSearch pairSearch,
std::unique_ptr< nbnxn_atomdata_t >  nbat,
const Nbnxm::KernelSetup &  kernelSetup,
std::unique_ptr< ExclusionChecker >  exclusionChecker,
NbnxmGpu *  gpu_nbv,
gmx_wallcycle *  wcycle 
)

Constructs an object from its components.

Parameters
[in]pairlistSetsThe pairlist sets, is consumed
[in]pairSearchThe pairsearch setup, is consumed
[in]nbatThe atom data, is consumed
[in]kernelSetupThe non-bonded kernel setup
[in]exclusionCheckerThe FEP exclusion checker, is consumed, can be nullptr
[in]gpu_nbvThe GPU non-bonded setup, ownership is transferred, can be nullptr
[in]wcyclePointer to wallcycle counters, can be nullptr
nonbonded_verlet_t::nonbonded_verlet_t ( std::unique_ptr< PairlistSets pairlistSets,
std::unique_ptr< PairSearch pairSearch,
std::unique_ptr< nbnxn_atomdata_t >  nbat,
const Nbnxm::KernelSetup &  kernelSetup,
NbnxmGpu *  gpu_nbv 
)

Constructs an object from its, minimal, components.

Parameters
[in]pairlistSetsThe pairlist sets, is consumed
[in]pairSearchThe pairsearch setup, is consumed
[in]nbatThe atom data, is consumed
[in]kernelSetupThe non-bonded kernel setup
[in]gpu_nbvThe GPU non-bonded setup, ownership is transferred, can be nullptr

Member Function Documentation

void nonbonded_verlet_t::atomdata_add_nbat_f_to_f ( gmx::AtomLocality  locality,
gmx::ArrayRef< gmx::RVec force 
)

Add the forces stored in nbat to f, zeros the forces in nbat.

Parameters
[in]localityLocal or non-local
[in,out]forceForce to be added to
void nonbonded_verlet_t::constructPairlist ( gmx::InteractionLocality  iLocality,
const gmx::ListOfLists< int > &  exclusions,
int64_t  step,
t_nrnb *  nrnb 
) const

Constructs the pairlist for the given locality.

When there are no non-self exclusions, exclusions can be empty. Otherwise the number of lists in exclusions should match the number of atoms when not using DD, or the total number of atoms in the i-zones when using DD.

Parameters
[in]iLocalityThe interaction locality: local or non-local
[in]exclusionsLists of exclusions for every atom.
[in]stepUsed to set the list creation step
[in,out]nrnbFlop accounting struct, can be nullptr
void nonbonded_verlet_t::convertCoordinates ( gmx::AtomLocality  locality,
gmx::ArrayRef< const gmx::RVec coordinates 
)

Convert the coordinates to NBNXM format for the given locality.

The API function for the transformation of the coordinates from one layout to another.

Parameters
[in]localityWhether coordinates for local or non-local atoms should be transformed.
[in]coordinatesCoordinates in plain rvec format to be transformed.
void nonbonded_verlet_t::convertCoordinatesGpu ( gmx::AtomLocality  locality,
DeviceBuffer< gmx::RVec d_x,
GpuEventSynchronizer *  xReadyOnDevice 
)

Convert the coordinates to NBNXM format on the GPU for the given locality.

The API function for the transformation of the coordinates from one layout to another in the GPU memory.

Parameters
[in]localityWhether coordinates for local or non-local atoms should be transformed.
[in]d_xGPU coordinates buffer in plain rvec format to be transformed.
[in]xReadyOnDeviceEvent synchronizer indicating that the coordinates are ready in the device memory.
int nonbonded_verlet_t::getNumAtoms ( gmx::AtomLocality  locality) const

Get the number of atoms for a given locality.

Parameters
[in]localityLocal or non-local
Returns
The number of atoms for given locality

The documentation for this struct was generated from the following files: