Gromacs
2021-beta3-UNCHECKED
|
#include <memory>
#include "gromacs/gpu_utils/devicebuffer_datatype.h"
#include "gromacs/math/vectypes.h"
#include "gromacs/mdtypes/locality.h"
#include "gromacs/utility/arrayref.h"
#include "gromacs/utility/enumerationhelpers.h"
#include "gromacs/utility/real.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 | |
class | gmx::ListOfLists< typename > |
A list of lists, optimized for performance. More... | |
class | gmx::Range< typename > |
Defines a range of integer numbers and accompanying operations. More... | |
struct | nonbonded_verlet_t |
Top-level non-bonded data structure for the Verlet-type cut-off scheme. More... | |
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, const t_inputrec *ir, const t_forcerec *fr, const t_commrec *cr, const gmx_hw_info_t &hardwareInfo, bool useGpuForNonbonded, const gmx::DeviceStreamManager *deviceStreamManager, 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 gridIndex, const rvec lowerCorner, const rvec upperCorner, const gmx::UpdateGroupsCog *updateGroupsCog, gmx::Range< int > atomRange, 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 | gridIndex, | ||
const rvec | lowerCorner, | ||
const rvec | upperCorner, | ||
const gmx::UpdateGroupsCog * | updateGroupsCog, | ||
gmx::Range< int > | atomRange, | ||
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 with indices wihtin atomRange
in x are put on the grid. When updateGroupsCog
!= nullptr, atoms are put on the grid based on the center of geometry of the group they belong to. Atoms or COGs of groups should be within the bounding box provided, this is checked in debug builds when not using update groups. The atom density is used to determine the grid size when gridIndex
= 0. When atomDensity
<= 0, the density is determined from atomEnd-atomStart and the bounding box corners. With domain decomposition, part of the atoms might have migrated, but have not been removed yet. This count is given by numAtomsMoved
. When move
[i] < 0 particle i has migrated and will not be put on the grid.
[in,out] | nb_verlet | The non-bonded object |
[in] | box | Box used for periodic distance calculations |
[in] | gridIndex | The index of the grid to spread to, always 0 except with test particle insertion |
[in] | lowerCorner | Atom groups to be gridded should have coordinates >= this corner |
[in] | upperCorner | Atom groups to be gridded should have coordinates <= this corner |
[in] | updateGroupsCog | Centers of geometry for update groups, pass nullptr when not using update groups |
[in] | atomRange | Range of atoms to grid |
[in] | atomDensity | An estimate of the atom density, used for peformance optimization and only with gridIndex = 0 |
[in] | atomInfo | Atom information flags |
[in] | x | Coordinates for atoms to grid |
[in] | numAtomsMoved | The number of atoms that will move to another domain, pass 0 without DD |
[in] | move | Move flags for atoms, pass nullptr without DD |
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.