Gromacs
2020-beta1
|
#include <memory>
#include "gromacs/gpu_utils/devicebuffer_datatype.h"
#include "gromacs/math/vectypes.h"
#include "gromacs/utility/arrayref.h"
#include "gromacs/utility/enumerationhelpers.h"
#include "gromacs/utility/real.h"
#include "locality.h"
#include "nbnxm_gpu.h"
This file contains the public interface of the nbnxm module that implements the NxM atom cluster non-bonded algorithm to efficiently compute pair forces.
Classes | |
struct | nonbonded_verlet_t |
Top-level non-bonded data structure for the Verlet-type cut-off scheme. More... | |
Enumerations | |
enum | BufferOpsUseGpu { True, False } |
Switch for whether to use GPU for buffer ops. | |
enum | KernelType : int { NotSet = 0, Cpu4x4_PlainC, Cpu4xN_Simd_4xN, Cpu4xN_Simd_2xNN, Gpu8x8x8, Cpu8x8x8_PlainC, Count } |
Nonbonded NxN kernel types: plain C, CPU SIMD, GPU, GPU emulation. | |
enum | EwaldExclusionType : int { NotSet = 0, Table, Analytical, DecidedByGpuModule } |
Ewald exclusion types. | |
enum | { enbvClearFNo, enbvClearFYes } |
Flag to tell the nonbonded kernels whether to clear the force output buffers. | |
Functions | |
const char * | Nbnxm::lookup_kernel_name (Nbnxm::KernelType kernelType) |
Return a string identifying the kernel type. More... | |
std::unique_ptr < nonbonded_verlet_t > | Nbnxm::init_nb_verlet (const gmx::MDLogger &mdlog, gmx_bool bFEP_NonBonded, const t_inputrec *ir, const t_forcerec *fr, const t_commrec *cr, const gmx_hw_info_t &hardwareInfo, const gmx_device_info_t *deviceInfo, const gmx_mtop_t *mtop, matrix box, gmx_wallcycle *wcycle) |
Creates an Nbnxm object. | |
void | nbnxn_put_on_grid (nonbonded_verlet_t *nb_verlet, const matrix box, int ddZone, const rvec lowerCorner, const rvec upperCorner, const gmx::UpdateGroupsCog *updateGroupsCog, int atomStart, int atomEnd, real atomDensity, gmx::ArrayRef< const int > atomInfo, gmx::ArrayRef< const gmx::RVec > x, int numAtomsMoved, const int *move) |
Put the atoms on the pair search grid. More... | |
void | nbnxn_put_on_grid_nonlocal (nonbonded_verlet_t *nb_verlet, const struct gmx_domdec_zones_t *zones, gmx::ArrayRef< const int > atomInfo, gmx::ArrayRef< const gmx::RVec > x) |
As nbnxn_put_on_grid, but for the non-local atoms. More... | |
void nbnxn_put_on_grid | ( | nonbonded_verlet_t * | nb_verlet, |
const matrix | box, | ||
int | ddZone, | ||
const rvec | lowerCorner, | ||
const rvec | upperCorner, | ||
const gmx::UpdateGroupsCog * | updateGroupsCog, | ||
int | atomStart, | ||
int | atomEnd, | ||
real | atomDensity, | ||
gmx::ArrayRef< const int > | atomInfo, | ||
gmx::ArrayRef< const gmx::RVec > | x, | ||
int | numAtomsMoved, | ||
const int * | move | ||
) |
Put the atoms on the pair search grid.
Only atoms atomStart to atomEnd in x are put on the grid. The atom_density is used to determine the grid size. When atomDensity<=0, the density is determined from atomEnd-atomStart and the corners. With domain decomposition part of the n particles might have migrated, but have not been removed yet. This count is given by nmoved. When move[i] < 0 particle i has migrated and will not be put on the grid. Without domain decomposition move will be NULL.
void nbnxn_put_on_grid_nonlocal | ( | nonbonded_verlet_t * | nb_verlet, |
const struct gmx_domdec_zones_t * | zones, | ||
gmx::ArrayRef< const int > | atomInfo, | ||
gmx::ArrayRef< const gmx::RVec > | x | ||
) |
As nbnxn_put_on_grid, but for the non-local atoms.
with domain decomposition. Should be called after calling nbnxn_search_put_on_grid for the local atoms / home zone.