Link GPU coordinate producer and consumer tasks
[alexxy/gromacs.git] / src / gromacs / nbnxm / nbnxm.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35
36 /*! \internal \file
37  * \brief
38  * Implements the Nbnxm class
39  *
40  * \author Berk Hess <hess@kth.se>
41  * \ingroup module_nbnxm
42  */
43
44 #include "gmxpre.h"
45
46 #include "nbnxm.h"
47
48 #include "gromacs/domdec/domdec_struct.h"
49 #include "gromacs/timing/wallcycle.h"
50
51 #include "atomdata.h"
52 #include "pairlistsets.h"
53 #include "pairsearch.h"
54
55 /*! \cond INTERNAL */
56
57 void nbnxn_put_on_grid(nonbonded_verlet_t             *nb_verlet,
58                        const matrix                    box,
59                        int                             ddZone,
60                        const rvec                      lowerCorner,
61                        const rvec                      upperCorner,
62                        const gmx::UpdateGroupsCog     *updateGroupsCog,
63                        int                             atomStart,
64                        int                             atomEnd,
65                        real                            atomDensity,
66                        gmx::ArrayRef<const int>        atomInfo,
67                        gmx::ArrayRef<const gmx::RVec>  x,
68                        int                             numAtomsMoved,
69                        const int                      *move)
70 {
71     nb_verlet->pairSearch_->putOnGrid(box, ddZone, lowerCorner, upperCorner,
72                                       updateGroupsCog, atomStart, atomEnd, atomDensity,
73                                       atomInfo, x, numAtomsMoved, move,
74                                       nb_verlet->nbat.get());
75 }
76
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)
82 {
83     for (int zone = 1; zone < zones->n; zone++)
84     {
85         rvec c0, c1;
86         for (int d = 0; d < DIM; d++)
87         {
88             c0[d] = zones->size[zone].bb_x0[d];
89             c1[d] = zones->size[zone].bb_x1[d];
90         }
91
92         nbnxn_put_on_grid(nbv, nullptr,
93                           zone, c0, c1,
94                           nullptr,
95                           zones->cg_range[zone],
96                           zones->cg_range[zone+1],
97                           -1,
98                           atomInfo,
99                           x,
100                           0, nullptr);
101     }
102 }
103
104 bool nonbonded_verlet_t::isDynamicPruningStepCpu(int64_t step) const
105 {
106     return pairlistSets_->isDynamicPruningStepCpu(step);
107 }
108
109 bool nonbonded_verlet_t::isDynamicPruningStepGpu(int64_t step) const
110 {
111     return pairlistSets_->isDynamicPruningStepGpu(step);
112 }
113
114 gmx::ArrayRef<const int> nonbonded_verlet_t::getLocalAtomOrder() const
115 {
116     /* Return the atom order for the home cell (index 0) */
117     const Nbnxm::Grid &grid       = pairSearch_->gridSet().grids()[0];
118
119     const int          numIndices = grid.atomIndexEnd() - grid.firstAtomInColumn(0);
120
121     return gmx::constArrayRefFromArray(pairSearch_->gridSet().atomIndices().data(), numIndices);
122 }
123
124 void nonbonded_verlet_t::setLocalAtomOrder()
125 {
126     pairSearch_->setLocalAtomOrder();
127 }
128
129 void nonbonded_verlet_t::setAtomProperties(const t_mdatoms          &mdatoms,
130                                            gmx::ArrayRef<const int>  atomInfo)
131 {
132     nbnxn_atomdata_set(nbat.get(), pairSearch_->gridSet(), &mdatoms, atomInfo.data());
133 }
134
135 void nonbonded_verlet_t::convertCoordinates(const Nbnxm::AtomLocality       locality,
136                                             const bool                      fillLocal,
137                                             gmx::ArrayRef<const gmx::RVec>  coordinates)
138 {
139     wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
140     wallcycle_sub_start(wcycle_, ewcsNB_X_BUF_OPS);
141
142     nbnxn_atomdata_copy_x_to_nbat_x(pairSearch_->gridSet(), locality, fillLocal,
143                                     as_rvec_array(coordinates.data()),
144                                     nbat.get());
145
146     wallcycle_sub_stop(wcycle_, ewcsNB_X_BUF_OPS);
147     wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
148 }
149
150 void nonbonded_verlet_t::convertCoordinatesGpu(const Nbnxm::AtomLocality        locality,
151                                                const bool                       fillLocal,
152                                                DeviceBuffer<float>              d_x,
153                                                GpuEventSynchronizer            *xReadyOnDevice)
154 {
155     wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
156     wallcycle_sub_start(wcycle_, ewcsNB_X_BUF_OPS);
157
158     nbnxn_atomdata_x_to_nbat_x_gpu(pairSearch_->gridSet(), locality, fillLocal,
159                                    gpu_nbv,
160                                    d_x,
161                                    xReadyOnDevice);
162
163     wallcycle_sub_stop(wcycle_, ewcsNB_X_BUF_OPS);
164     wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
165 }
166
167 gmx::ArrayRef<const int> nonbonded_verlet_t::getGridIndices() const
168 {
169     return pairSearch_->gridSet().cells();
170 }
171
172 void
173 nonbonded_verlet_t::atomdata_add_nbat_f_to_f(const Nbnxm::AtomLocality           locality,
174                                              gmx::ArrayRef<gmx::RVec>            force)
175 {
176
177     /* Skip the reduction if there was no short-range GPU work to do
178      * (either NB or both NB and bonded work). */
179     if (!pairlistIsSimple() && !haveGpuShortRangeWork(locality))
180     {
181         return;
182     }
183
184     wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
185     wallcycle_sub_start(wcycle_, ewcsNB_F_BUF_OPS);
186
187     reduceForces(nbat.get(), locality, pairSearch_->gridSet(), as_rvec_array(force.data()));
188
189     wallcycle_sub_stop(wcycle_, ewcsNB_F_BUF_OPS);
190     wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
191 }
192
193 void
194 nonbonded_verlet_t::atomdata_add_nbat_f_to_f_gpu(const Nbnxm::AtomLocality                   locality,
195                                                  DeviceBuffer<float>                         totalForcesDevice,
196                                                  void                                       *forcesPmeDevice,
197                                                  gmx::ArrayRef<GpuEventSynchronizer* const>  dependencyList,
198                                                  bool                                        useGpuFPmeReduction,
199                                                  bool                                        accumulateForce)
200 {
201
202     GMX_ASSERT((useGpuFPmeReduction == (forcesPmeDevice != nullptr)),
203                "GPU PME force reduction is only valid when a non-null GPU PME force pointer is available");
204
205     /* Skip the reduction if there was no short-range GPU work to do
206      * (either NB or both NB and bonded work). */
207     if (!pairlistIsSimple() && !haveGpuShortRangeWork(locality))
208     {
209         return;
210     }
211
212     wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
213     wallcycle_sub_start(wcycle_, ewcsNB_F_BUF_OPS);
214
215     reduceForcesGpu(locality, totalForcesDevice, pairSearch_->gridSet(), forcesPmeDevice, dependencyList, gpu_nbv, useGpuFPmeReduction, accumulateForce);
216
217     wallcycle_sub_stop(wcycle_, ewcsNB_F_BUF_OPS);
218     wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
219 }
220
221 void
222 nonbonded_verlet_t::atomdata_init_add_nbat_f_to_f_gpu()
223 {
224
225     wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
226     wallcycle_sub_start(wcycle_, ewcsNB_F_BUF_OPS);
227
228     const Nbnxm::GridSet      &gridSet = pairSearch_->gridSet();
229
230     Nbnxm::nbnxn_gpu_init_add_nbat_f_to_f(gridSet.cells().data(),
231                                           gpu_nbv,
232                                           gridSet.numRealAtomsTotal());
233
234     wallcycle_sub_stop(wcycle_, ewcsNB_F_BUF_OPS);
235     wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
236 }
237
238 real nonbonded_verlet_t::pairlistInnerRadius() const
239 {
240     return pairlistSets_->params().rlistInner;
241 }
242
243 real nonbonded_verlet_t::pairlistOuterRadius() const
244 {
245     return pairlistSets_->params().rlistOuter;
246 }
247
248 void nonbonded_verlet_t::changePairlistRadii(real rlistOuter,
249                                              real rlistInner)
250 {
251     pairlistSets_->changePairlistRadii(rlistOuter, rlistInner);
252 }
253
254 void
255 nonbonded_verlet_t::atomdata_init_copy_x_to_nbat_x_gpu()
256 {
257     Nbnxm::nbnxn_gpu_init_x_to_nbat_x(pairSearch_->gridSet(), gpu_nbv);
258 }
259
260 void nonbonded_verlet_t::insertNonlocalGpuDependency(const Nbnxm::InteractionLocality interactionLocality)
261 {
262     Nbnxm::nbnxnInsertNonlocalGpuDependency(gpu_nbv, interactionLocality);
263 }
264
265 void* nonbonded_verlet_t::get_x_on_device_event()
266 {
267     return Nbnxm::nbnxn_get_x_on_device_event(gpu_nbv);
268 }
269
270 void nonbonded_verlet_t::wait_nonlocal_x_copy_D2H_done()
271 {
272     Nbnxm::nbnxn_wait_nonlocal_x_copy_D2H_done(gpu_nbv);
273 }
274
275 void nonbonded_verlet_t::stream_local_wait_for_nonlocal()
276 {
277     Nbnxm::nbnxn_stream_local_wait_for_nonlocal(gpu_nbv);
278 }
279
280 /*! \endcond */