2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019,2020,2021, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
38 * Implements the Nbnxm class
40 * \author Berk Hess <hess@kth.se>
41 * \ingroup module_nbnxm
48 #include "gromacs/domdec/domdec_struct.h"
49 #include "gromacs/nbnxm/atomdata.h"
50 #include "gromacs/timing/wallcycle.h"
52 #include "nbnxm_gpu.h"
53 #include "pairlistsets.h"
54 #include "pairsearch.h"
58 void nbnxn_put_on_grid(nonbonded_verlet_t* nb_verlet,
61 const rvec lowerCorner,
62 const rvec upperCorner,
63 const gmx::UpdateGroupsCog* updateGroupsCog,
64 gmx::Range<int> atomRange,
66 gmx::ArrayRef<const int> atomInfo,
67 gmx::ArrayRef<const gmx::RVec> x,
71 nb_verlet->pairSearch_->putOnGrid(box,
82 nb_verlet->nbat.get());
85 /* Calls nbnxn_put_on_grid for all non-local domains */
86 void nbnxn_put_on_grid_nonlocal(nonbonded_verlet_t* nbv,
87 const struct gmx_domdec_zones_t* zones,
88 gmx::ArrayRef<const int> atomInfo,
89 gmx::ArrayRef<const gmx::RVec> x)
91 for (int zone = 1; zone < zones->n; zone++)
94 for (int d = 0; d < DIM; d++)
96 c0[d] = zones->size[zone].bb_x0[d];
97 c1[d] = zones->size[zone].bb_x1[d];
100 nbnxn_put_on_grid(nbv,
106 { zones->cg_range[zone], zones->cg_range[zone + 1] },
115 bool nonbonded_verlet_t::isDynamicPruningStepCpu(int64_t step) const
117 return pairlistSets_->isDynamicPruningStepCpu(step);
120 bool nonbonded_verlet_t::isDynamicPruningStepGpu(int64_t step) const
122 return pairlistSets_->isDynamicPruningStepGpu(step);
125 gmx::ArrayRef<const int> nonbonded_verlet_t::getLocalAtomOrder() const
127 /* Return the atom order for the home cell (index 0) */
128 const Nbnxm::Grid& grid = pairSearch_->gridSet().grids()[0];
130 const int numIndices = grid.atomIndexEnd() - grid.firstAtomInColumn(0);
132 return gmx::constArrayRefFromArray(pairSearch_->gridSet().atomIndices().data(), numIndices);
135 void nonbonded_verlet_t::setLocalAtomOrder() const
137 pairSearch_->setLocalAtomOrder();
140 void nonbonded_verlet_t::setAtomProperties(gmx::ArrayRef<const int> atomTypes,
141 gmx::ArrayRef<const real> atomCharges,
142 gmx::ArrayRef<const int> atomInfo) const
144 nbnxn_atomdata_set(nbat.get(), pairSearch_->gridSet(), atomTypes, atomCharges, atomInfo);
147 void nonbonded_verlet_t::convertCoordinates(const gmx::AtomLocality locality,
148 const bool fillLocal,
149 gmx::ArrayRef<const gmx::RVec> coordinates)
151 wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
152 wallcycle_sub_start(wcycle_, ewcsNB_X_BUF_OPS);
154 nbnxn_atomdata_copy_x_to_nbat_x(
155 pairSearch_->gridSet(), locality, fillLocal, as_rvec_array(coordinates.data()), nbat.get());
157 wallcycle_sub_stop(wcycle_, ewcsNB_X_BUF_OPS);
158 wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
161 void nonbonded_verlet_t::convertCoordinatesGpu(const gmx::AtomLocality locality,
162 const bool fillLocal,
163 DeviceBuffer<gmx::RVec> d_x,
164 GpuEventSynchronizer* xReadyOnDevice)
166 wallcycle_start(wcycle_, ewcLAUNCH_GPU);
167 wallcycle_sub_start(wcycle_, ewcsLAUNCH_GPU_NB_X_BUF_OPS);
169 nbnxn_atomdata_x_to_nbat_x_gpu(pairSearch_->gridSet(), locality, fillLocal, gpu_nbv, d_x, xReadyOnDevice);
171 wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_NB_X_BUF_OPS);
172 wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
175 gmx::ArrayRef<const int> nonbonded_verlet_t::getGridIndices() const
177 return pairSearch_->gridSet().cells();
180 void nonbonded_verlet_t::atomdata_add_nbat_f_to_f(const gmx::AtomLocality locality,
181 gmx::ArrayRef<gmx::RVec> force)
184 /* Skip the reduction if there was no short-range GPU work to do
185 * (either NB or both NB and bonded work). */
186 if (!pairlistIsSimple() && !Nbnxm::haveGpuShortRangeWork(gpu_nbv, locality))
191 wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
192 wallcycle_sub_start(wcycle_, ewcsNB_F_BUF_OPS);
194 reduceForces(nbat.get(), locality, pairSearch_->gridSet(), as_rvec_array(force.data()));
196 wallcycle_sub_stop(wcycle_, ewcsNB_F_BUF_OPS);
197 wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
200 int nonbonded_verlet_t::getNumAtoms(const gmx::AtomLocality locality) const
205 case gmx::AtomLocality::All: numAtoms = pairSearch_->gridSet().numRealAtomsTotal(); break;
206 case gmx::AtomLocality::Local: numAtoms = pairSearch_->gridSet().numRealAtomsLocal(); break;
207 case gmx::AtomLocality::NonLocal:
208 numAtoms = pairSearch_->gridSet().numRealAtomsTotal()
209 - pairSearch_->gridSet().numRealAtomsLocal();
211 case gmx::AtomLocality::Count:
212 GMX_ASSERT(false, "Count is invalid locality specifier");
218 void* nonbonded_verlet_t::getGpuForces() const
220 return Nbnxm::getGpuForces(gpu_nbv);
223 real nonbonded_verlet_t::pairlistInnerRadius() const
225 return pairlistSets_->params().rlistInner;
228 real nonbonded_verlet_t::pairlistOuterRadius() const
230 return pairlistSets_->params().rlistOuter;
233 void nonbonded_verlet_t::changePairlistRadii(real rlistOuter, real rlistInner) const
235 pairlistSets_->changePairlistRadii(rlistOuter, rlistInner);
238 void nonbonded_verlet_t::setupGpuShortRangeWork(const gmx::GpuBonded* gpuBonded,
239 const gmx::InteractionLocality iLocality) const
241 if (useGpu() && !emulateGpu())
243 Nbnxm::setupGpuShortRangeWork(gpu_nbv, gpuBonded, iLocality);
247 void nonbonded_verlet_t::atomdata_init_copy_x_to_nbat_x_gpu() const
249 Nbnxm::nbnxn_gpu_init_x_to_nbat_x(pairSearch_->gridSet(), gpu_nbv);