Gromacs
2025.0-dev-20241014-f673b97
|
#include <gromacs/nbnxm/atomdata.h>
Struct that stores atom related data for the nbnxn module.
Note: performance would improve slightly when all std::vector containers in this struct would not initialize during resize().
Classes | |
struct | Params |
The actual atom data parameter values. More... | |
struct | SimdMasks |
Diagonal and topology exclusion helper data for all SIMD kernels. More... | |
Public Member Functions | |
nbnxn_atomdata_t (PinningPolicy pinningPolicy) | |
nbnxn_atomdata_t (PinningPolicy pinningPolicy, const MDLogger &mdlog, NbnxmKernelType kernelType, const std::optional< LJCombinationRule > &ljCombinationRule, LJCombinationRule pmeLJCombinationRule, ArrayRef< const real > nbfp, bool addFillerAtomType, int numEnergyGroups, int numOutputBuffers) | |
Constructor. More... | |
const Params & | params () const |
Returns a const reference to the parameters. | |
Params & | paramsDeprecated () |
Returns a non-const reference to the parameters. | |
int | numAtoms () const |
Returns the current total number of atoms stored. | |
int | numLocalAtoms () const |
Returns the number of local atoms. | |
ArrayRef< const real > | x () const |
Return the coordinate buffer, and q with xFormat==nbatXYZQ. | |
ArrayRef< real > | x () |
Return the coordinate buffer, and q with xFormat==nbatXYZQ. | |
const SimdMasks & | simdMasks () const |
Masks for handling exclusions in the SIMD kernels. | |
void | resizeCoordinateBuffer (int numAtoms, int domainDecompositionZone) |
Resizes the coordinate buffer and sets the number of atoms. More... | |
void | resizeForceBuffers () |
Resizes the force buffers for the current number of atoms. | |
nbnxn_atomdata_output_t & | outputBuffer (int thread) |
Returns the output buffer for the given thread . | |
const nbnxn_atomdata_output_t & | outputBuffer (int thread) const |
Returns the output buffer for the given thread . | |
ArrayRef< const nbnxn_atomdata_output_t > | outputBuffers () const |
Returns the list of output buffers. | |
bool | useBufferFlags () const |
Returns whether buffer flags are used. | |
std::vector< gmx_bitmask_t > & | bufferFlags () |
Returns the vector of buffer flags. | |
void | reduceForces (AtomLocality locality, const GridSet &gridSet, ArrayRef< RVec > totalForce) |
Add the computed forces to f , an internal reduction might be performed as well. More... | |
void | clearForceBuffer (int outputIndex) |
Clears the force buffer. More... | |
Public Attributes | |
int | XFormat |
The format of x (and q), enum. | |
int | FFormat |
The format of f, enum. | |
bool | bDynamicBox |
Do we need to update shift_vec every step? | |
HostVector< RVec > | shift_vec |
Shift vectors, copied from t_forcerec. | |
int | xstride |
stride for a coordinate in x (usually 3 or 4) | |
int | fstride |
stride for a coordinate in f (usually 3 or 4) | |
gmx::nbnxn_atomdata_t::nbnxn_atomdata_t | ( | PinningPolicy | pinningPolicy, |
const MDLogger & | mdlog, | ||
NbnxmKernelType | kernelType, | ||
const std::optional< LJCombinationRule > & | ljCombinationRule, | ||
LJCombinationRule | pmeLJCombinationRule, | ||
ArrayRef< const real > | nbfp, | ||
bool | addFillerAtomType, | ||
int | numEnergyGroups, | ||
int | numOutputBuffers | ||
) |
Constructor.
This class only stores one, or no, LJ combination rule parameter list. With LJ-PME the rule for the LJ PME-grid must match the PME grid combination rule and there can be no combination rule for the LJ pair parameters. Without LJ-PME the combination rule should match the combination rule for the LJ parameters. An (release) assertion failure will occur when these conditions are not met.
The non-bonded force parameter matrix nbfp
is a serialized version of a square parameter matrix with a pair of parameters 6*C6, 12*C12 for every atom type pair.
[in] | pinningPolicy | Sets the pinning policy for all data that might be transferred to a GPU |
[in] | mdlog | The logger |
[in] | kernelType | Nonbonded NxN kernel type |
[in] | ljCombinationRule | The LJ combination rule parameters to generate, empty is detect from the LJ parameters |
[in] | pmeLJCombinationRule | The LJ combination rule parameters to generate for the LJ PME-grid part |
[in] | nbfp | Non-bonded force parameter matrix |
[in] | addFillerAtomType | When true, add a filler atom type, when false, nbfp should have atom-type with index numTypes-1 with all parameters zero so that that row and column have only zero values. |
[in] | numEnergyGroups | Number of energy groups |
[in] | numOutputBuffers | Number of output data structures |
void gmx::nbnxn_atomdata_t::clearForceBuffer | ( | int | outputIndex | ) |
Clears the force buffer.
Either the whole buffer is cleared or only the parts used by thread/task outputIndex
when useBufferFlags() returns true.
[in] | outputIndex | The index of the output object to clear |
void gmx::nbnxn_atomdata_t::reduceForces | ( | AtomLocality | locality, |
const GridSet & | gridSet, | ||
ArrayRef< RVec > | totalForce | ||
) |
Add the computed forces to f
, an internal reduction might be performed as well.
[in] | locality | If the reduction should be performed on local or non-local atoms. |
[in] | gridSet | The grids data. |
[out] | totalForce | Buffer to accumulate resulting force |
void gmx::nbnxn_atomdata_t::resizeCoordinateBuffer | ( | int | numAtoms, |
int | domainDecompositionZone | ||
) |
Resizes the coordinate buffer and sets the number of atoms.
numAtoms | The new number of atoms |
domainDecompositionZone | The domain decomposition zone index |