Update instructions in containers.rst
[gromacs.git] / src / gromacs / nbnxm / nbnxm.cpp
blob76e2eead9f07d18144c568948f3b7b5c34b95667
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019,2020, 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.
36 /*! \internal \file
37 * \brief
38 * Implements the Nbnxm class
40 * \author Berk Hess <hess@kth.se>
41 * \ingroup module_nbnxm
44 #include "gmxpre.h"
46 #include "nbnxm.h"
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"
56 /*! \cond INTERNAL */
58 void nbnxn_put_on_grid(nonbonded_verlet_t* nb_verlet,
59 const matrix box,
60 int gridIndex,
61 const rvec lowerCorner,
62 const rvec upperCorner,
63 const gmx::UpdateGroupsCog* updateGroupsCog,
64 gmx::Range<int> atomRange,
65 real atomDensity,
66 gmx::ArrayRef<const int> atomInfo,
67 gmx::ArrayRef<const gmx::RVec> x,
68 int numAtomsMoved,
69 const int* move)
71 nb_verlet->pairSearch_->putOnGrid(box, gridIndex, lowerCorner, upperCorner, updateGroupsCog,
72 atomRange, atomDensity, atomInfo, x, numAtomsMoved, move,
73 nb_verlet->nbat.get());
76 /* Calls nbnxn_put_on_grid for all non-local domains */
77 void nbnxn_put_on_grid_nonlocal(nonbonded_verlet_t* nbv,
78 const struct gmx_domdec_zones_t* zones,
79 gmx::ArrayRef<const int> atomInfo,
80 gmx::ArrayRef<const gmx::RVec> x)
82 for (int zone = 1; zone < zones->n; zone++)
84 rvec c0, c1;
85 for (int d = 0; d < DIM; d++)
87 c0[d] = zones->size[zone].bb_x0[d];
88 c1[d] = zones->size[zone].bb_x1[d];
91 nbnxn_put_on_grid(nbv, nullptr, zone, c0, c1, nullptr,
92 { zones->cg_range[zone], zones->cg_range[zone + 1] }, -1, atomInfo, x, 0,
93 nullptr);
97 bool nonbonded_verlet_t::isDynamicPruningStepCpu(int64_t step) const
99 return pairlistSets_->isDynamicPruningStepCpu(step);
102 bool nonbonded_verlet_t::isDynamicPruningStepGpu(int64_t step) const
104 return pairlistSets_->isDynamicPruningStepGpu(step);
107 gmx::ArrayRef<const int> nonbonded_verlet_t::getLocalAtomOrder() const
109 /* Return the atom order for the home cell (index 0) */
110 const Nbnxm::Grid& grid = pairSearch_->gridSet().grids()[0];
112 const int numIndices = grid.atomIndexEnd() - grid.firstAtomInColumn(0);
114 return gmx::constArrayRefFromArray(pairSearch_->gridSet().atomIndices().data(), numIndices);
117 void nonbonded_verlet_t::setLocalAtomOrder()
119 pairSearch_->setLocalAtomOrder();
122 void nonbonded_verlet_t::setAtomProperties(gmx::ArrayRef<const int> atomTypes,
123 gmx::ArrayRef<const real> atomCharges,
124 gmx::ArrayRef<const int> atomInfo)
126 nbnxn_atomdata_set(nbat.get(), pairSearch_->gridSet(), atomTypes, atomCharges, atomInfo);
129 void nonbonded_verlet_t::convertCoordinates(const gmx::AtomLocality locality,
130 const bool fillLocal,
131 gmx::ArrayRef<const gmx::RVec> coordinates)
133 wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
134 wallcycle_sub_start(wcycle_, ewcsNB_X_BUF_OPS);
136 nbnxn_atomdata_copy_x_to_nbat_x(pairSearch_->gridSet(), locality, fillLocal,
137 as_rvec_array(coordinates.data()), nbat.get());
139 wallcycle_sub_stop(wcycle_, ewcsNB_X_BUF_OPS);
140 wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
143 void nonbonded_verlet_t::convertCoordinatesGpu(const gmx::AtomLocality locality,
144 const bool fillLocal,
145 DeviceBuffer<gmx::RVec> d_x,
146 GpuEventSynchronizer* xReadyOnDevice)
148 wallcycle_start(wcycle_, ewcLAUNCH_GPU);
149 wallcycle_sub_start(wcycle_, ewcsLAUNCH_GPU_NB_X_BUF_OPS);
151 nbnxn_atomdata_x_to_nbat_x_gpu(pairSearch_->gridSet(), locality, fillLocal, gpu_nbv, d_x, xReadyOnDevice);
153 wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_NB_X_BUF_OPS);
154 wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
157 gmx::ArrayRef<const int> nonbonded_verlet_t::getGridIndices() const
159 return pairSearch_->gridSet().cells();
162 void nonbonded_verlet_t::atomdata_add_nbat_f_to_f(const gmx::AtomLocality locality,
163 gmx::ArrayRef<gmx::RVec> force)
166 /* Skip the reduction if there was no short-range GPU work to do
167 * (either NB or both NB and bonded work). */
168 if (!pairlistIsSimple() && !Nbnxm::haveGpuShortRangeWork(gpu_nbv, locality))
170 return;
173 wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
174 wallcycle_sub_start(wcycle_, ewcsNB_F_BUF_OPS);
176 reduceForces(nbat.get(), locality, pairSearch_->gridSet(), as_rvec_array(force.data()));
178 wallcycle_sub_stop(wcycle_, ewcsNB_F_BUF_OPS);
179 wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
182 int nonbonded_verlet_t::getNumAtoms(const gmx::AtomLocality locality)
184 int numAtoms = 0;
185 switch (locality)
187 case gmx::AtomLocality::All: numAtoms = pairSearch_->gridSet().numRealAtomsTotal(); break;
188 case gmx::AtomLocality::Local: numAtoms = pairSearch_->gridSet().numRealAtomsLocal(); break;
189 case gmx::AtomLocality::NonLocal:
190 numAtoms = pairSearch_->gridSet().numRealAtomsTotal()
191 - pairSearch_->gridSet().numRealAtomsLocal();
192 break;
193 case gmx::AtomLocality::Count:
194 GMX_ASSERT(false, "Count is invalid locality specifier");
195 break;
197 return numAtoms;
200 void* nonbonded_verlet_t::getGpuForces()
202 return Nbnxm::getGpuForces(gpu_nbv);
205 real nonbonded_verlet_t::pairlistInnerRadius() const
207 return pairlistSets_->params().rlistInner;
210 real nonbonded_verlet_t::pairlistOuterRadius() const
212 return pairlistSets_->params().rlistOuter;
215 void nonbonded_verlet_t::changePairlistRadii(real rlistOuter, real rlistInner)
217 pairlistSets_->changePairlistRadii(rlistOuter, rlistInner);
220 void nonbonded_verlet_t::setupGpuShortRangeWork(const gmx::GpuBonded* gpuBonded,
221 const gmx::InteractionLocality iLocality)
223 if (useGpu() && !emulateGpu())
225 Nbnxm::setupGpuShortRangeWork(gpu_nbv, gpuBonded, iLocality);
229 void nonbonded_verlet_t::atomdata_init_copy_x_to_nbat_x_gpu()
231 Nbnxm::nbnxn_gpu_init_x_to_nbat_x(pairSearch_->gridSet(), gpu_nbv);
234 void nonbonded_verlet_t::insertNonlocalGpuDependency(const gmx::InteractionLocality interactionLocality)
236 Nbnxm::nbnxnInsertNonlocalGpuDependency(gpu_nbv, interactionLocality);
239 /*! \endcond */