2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019, 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 "pairlistsets.h"
53 #include "pairsearch.h"
57 void nbnxn_put_on_grid(nonbonded_verlet_t *nb_verlet,
60 const rvec lowerCorner,
61 const rvec upperCorner,
62 const gmx::UpdateGroupsCog *updateGroupsCog,
66 gmx::ArrayRef<const int> atomInfo,
67 gmx::ArrayRef<const gmx::RVec> x,
71 nb_verlet->pairSearch_->putOnGrid(box, ddZone, lowerCorner, upperCorner,
72 updateGroupsCog, atomStart, atomEnd, atomDensity,
73 atomInfo, x, numAtomsMoved, move,
74 nb_verlet->nbat.get());
77 /* Calls nbnxn_put_on_grid for all non-local domains */
78 void nbnxn_put_on_grid_nonlocal(nonbonded_verlet_t *nbv,
79 const struct gmx_domdec_zones_t *zones,
80 gmx::ArrayRef<const int> atomInfo,
81 gmx::ArrayRef<const gmx::RVec> x)
83 for (int zone = 1; zone < zones->n; zone++)
86 for (int d = 0; d < DIM; d++)
88 c0[d] = zones->size[zone].bb_x0[d];
89 c1[d] = zones->size[zone].bb_x1[d];
92 nbnxn_put_on_grid(nbv, nullptr,
95 zones->cg_range[zone],
96 zones->cg_range[zone+1],
104 bool nonbonded_verlet_t::isDynamicPruningStepCpu(int64_t step) const
106 return pairlistSets_->isDynamicPruningStepCpu(step);
109 bool nonbonded_verlet_t::isDynamicPruningStepGpu(int64_t step) const
111 return pairlistSets_->isDynamicPruningStepGpu(step);
114 gmx::ArrayRef<const int> nonbonded_verlet_t::getLocalAtomOrder() const
116 /* Return the atom order for the home cell (index 0) */
117 const Nbnxm::Grid &grid = pairSearch_->gridSet().grids()[0];
119 const int numIndices = grid.atomIndexEnd() - grid.firstAtomInColumn(0);
121 return gmx::constArrayRefFromArray(pairSearch_->gridSet().atomIndices().data(), numIndices);
124 void nonbonded_verlet_t::setLocalAtomOrder()
126 pairSearch_->setLocalAtomOrder();
129 void nonbonded_verlet_t::setAtomProperties(const t_mdatoms &mdatoms,
130 gmx::ArrayRef<const int> atomInfo)
132 nbnxn_atomdata_set(nbat.get(), pairSearch_->gridSet(), &mdatoms, atomInfo.data());
135 void nonbonded_verlet_t::convertCoordinates(const Nbnxm::AtomLocality locality,
136 const bool fillLocal,
137 gmx::ArrayRef<const gmx::RVec> coordinates)
139 wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
140 wallcycle_sub_start(wcycle_, ewcsNB_X_BUF_OPS);
142 nbnxn_atomdata_copy_x_to_nbat_x(pairSearch_->gridSet(), locality, fillLocal,
143 as_rvec_array(coordinates.data()),
146 wallcycle_sub_stop(wcycle_, ewcsNB_X_BUF_OPS);
147 wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
151 void nonbonded_verlet_t::copyCoordinatesToGpu(const Nbnxm::AtomLocality locality,
152 const bool fillLocal,
153 gmx::ArrayRef<const gmx::RVec> coordinatesHost)
155 wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
156 wallcycle_sub_start(wcycle_, ewcsNB_X_BUF_OPS);
158 nbnxn_atomdata_copy_x_to_gpu(pairSearch_->gridSet(), locality, fillLocal,
160 as_rvec_array(coordinatesHost.data()));
162 wallcycle_sub_stop(wcycle_, ewcsNB_X_BUF_OPS);
163 wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
166 DeviceBuffer<float> nonbonded_verlet_t::getDeviceCoordinates()
168 return nbnxn_atomdata_get_x_gpu(gpu_nbv);
171 void nonbonded_verlet_t::convertCoordinatesGpu(const Nbnxm::AtomLocality locality,
172 const bool fillLocal,
173 DeviceBuffer<float> coordinatesDevice)
175 wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
176 wallcycle_sub_start(wcycle_, ewcsNB_X_BUF_OPS);
178 nbnxn_atomdata_x_to_nbat_x_gpu(pairSearch_->gridSet(), locality, fillLocal,
182 wallcycle_sub_stop(wcycle_, ewcsNB_X_BUF_OPS);
183 wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
186 gmx::ArrayRef<const int> nonbonded_verlet_t::getGridIndices() const
188 return pairSearch_->gridSet().cells();
192 nonbonded_verlet_t::atomdata_add_nbat_f_to_f(const Nbnxm::AtomLocality locality,
193 gmx::ArrayRef<gmx::RVec> force)
196 /* Skip the reduction if there was no short-range GPU work to do
197 * (either NB or both NB and bonded work). */
198 if (!pairlistIsSimple() && !haveGpuShortRangeWork(locality))
203 wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
204 wallcycle_sub_start(wcycle_, ewcsNB_F_BUF_OPS);
206 reduceForces(nbat.get(), locality, pairSearch_->gridSet(), as_rvec_array(force.data()));
208 wallcycle_sub_stop(wcycle_, ewcsNB_F_BUF_OPS);
209 wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
213 nonbonded_verlet_t::atomdata_add_nbat_f_to_f_gpu(const Nbnxm::AtomLocality locality,
214 DeviceBuffer<float> totalForcesDevice,
215 void *forcesPmeDevice,
216 GpuEventSynchronizer *pmeForcesReady,
217 bool useGpuFPmeReduction,
218 bool accumulateForce)
221 GMX_ASSERT((useGpuFPmeReduction == (forcesPmeDevice != nullptr)),
222 "GPU PME force reduction is only valid when a non-null GPU PME force pointer is available");
224 /* Skip the reduction if there was no short-range GPU work to do
225 * (either NB or both NB and bonded work). */
226 if (!pairlistIsSimple() && !haveGpuShortRangeWork(locality))
231 wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
232 wallcycle_sub_start(wcycle_, ewcsNB_F_BUF_OPS);
234 reduceForcesGpu(locality, totalForcesDevice, pairSearch_->gridSet(), forcesPmeDevice, pmeForcesReady, gpu_nbv, useGpuFPmeReduction, accumulateForce);
236 wallcycle_sub_stop(wcycle_, ewcsNB_F_BUF_OPS);
237 wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
241 nonbonded_verlet_t::atomdata_init_add_nbat_f_to_f_gpu()
244 wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
245 wallcycle_sub_start(wcycle_, ewcsNB_F_BUF_OPS);
247 const Nbnxm::GridSet &gridSet = pairSearch_->gridSet();
249 Nbnxm::nbnxn_gpu_init_add_nbat_f_to_f(gridSet.cells().data(),
251 gridSet.numRealAtomsTotal());
253 wallcycle_sub_stop(wcycle_, ewcsNB_F_BUF_OPS);
254 wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
257 DeviceBuffer<float> nonbonded_verlet_t::getDeviceForces()
259 return nbnxn_atomdata_get_f_gpu(gpu_nbv);
262 real nonbonded_verlet_t::pairlistInnerRadius() const
264 return pairlistSets_->params().rlistInner;
267 real nonbonded_verlet_t::pairlistOuterRadius() const
269 return pairlistSets_->params().rlistOuter;
272 void nonbonded_verlet_t::changePairlistRadii(real rlistOuter,
275 pairlistSets_->changePairlistRadii(rlistOuter, rlistInner);
279 nonbonded_verlet_t::atomdata_init_copy_x_to_nbat_x_gpu()
281 Nbnxm::nbnxn_gpu_init_x_to_nbat_x(pairSearch_->gridSet(), gpu_nbv);
284 void nonbonded_verlet_t::insertNonlocalGpuDependency(const Nbnxm::InteractionLocality interactionLocality)
286 Nbnxm::nbnxnInsertNonlocalGpuDependency(gpu_nbv, interactionLocality);
289 void nonbonded_verlet_t::launch_copy_f_to_gpu(rvec *f, const Nbnxm::AtomLocality locality)
291 nbnxn_launch_copy_f_to_gpu(locality,
292 pairSearch_->gridSet(),
297 void nonbonded_verlet_t::launch_copy_f_from_gpu(rvec *f, const Nbnxm::AtomLocality locality)
299 nbnxn_launch_copy_f_from_gpu(locality,
300 pairSearch_->gridSet(),
305 void nonbonded_verlet_t::wait_for_gpu_force_reduction(const Nbnxm::AtomLocality locality)
307 nbnxn_wait_for_gpu_force_reduction(locality, gpu_nbv);