Decouple coordinates buffer management from buffer ops in NBNXM
[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/nbnxm/atomdata.h"
50 #include "gromacs/timing/wallcycle.h"
51
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
151 void nonbonded_verlet_t::copyCoordinatesToGpu(const Nbnxm::AtomLocality       locality,
152                                               const bool                      fillLocal,
153                                               gmx::ArrayRef<const gmx::RVec>  coordinatesHost)
154 {
155     wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
156     wallcycle_sub_start(wcycle_, ewcsNB_X_BUF_OPS);
157
158     nbnxn_atomdata_copy_x_to_gpu(pairSearch_->gridSet(), locality, fillLocal,
159                                  nbat.get(), gpu_nbv,
160                                  as_rvec_array(coordinatesHost.data()));
161
162     wallcycle_sub_stop(wcycle_, ewcsNB_X_BUF_OPS);
163     wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
164 }
165
166 DeviceBuffer<float> nonbonded_verlet_t::getDeviceCoordinates()
167 {
168     return nbnxn_atomdata_get_x_gpu(gpu_nbv);
169 }
170
171 void nonbonded_verlet_t::convertCoordinatesGpu(const Nbnxm::AtomLocality       locality,
172                                                const bool                      fillLocal,
173                                                DeviceBuffer<float>             coordinatesDevice)
174 {
175     wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
176     wallcycle_sub_start(wcycle_, ewcsNB_X_BUF_OPS);
177
178     nbnxn_atomdata_x_to_nbat_x_gpu(pairSearch_->gridSet(), locality, fillLocal,
179                                    gpu_nbv,
180                                    coordinatesDevice);
181
182     wallcycle_sub_stop(wcycle_, ewcsNB_X_BUF_OPS);
183     wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
184 }
185
186 gmx::ArrayRef<const int> nonbonded_verlet_t::getGridIndices() const
187 {
188     return pairSearch_->gridSet().cells();
189 }
190
191 void
192 nonbonded_verlet_t::atomdata_add_nbat_f_to_f(const Nbnxm::AtomLocality           locality,
193                                              gmx::ArrayRef<gmx::RVec>            force)
194 {
195
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))
199     {
200         return;
201     }
202
203     wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
204     wallcycle_sub_start(wcycle_, ewcsNB_F_BUF_OPS);
205
206     reduceForces<false>(nbat.get(), locality, pairSearch_->gridSet(), as_rvec_array(force.data()), nullptr, nullptr, gpu_nbv, false, false);
207
208     wallcycle_sub_stop(wcycle_, ewcsNB_F_BUF_OPS);
209     wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
210 }
211
212 void
213 nonbonded_verlet_t::atomdata_add_nbat_f_to_f(const Nbnxm::AtomLocality           locality,
214                                              gmx::ArrayRef<gmx::RVec>            force,
215                                              void                               *fPmeDeviceBuffer,
216                                              GpuEventSynchronizer               *pmeForcesReady,
217                                              BufferOpsUseGpu                     useGpu,
218                                              bool                                useGpuFPmeReduction,
219                                              bool                                accumulateForce)
220 {
221
222     GMX_ASSERT(!((useGpu == BufferOpsUseGpu::False) && accumulateForce),
223                "Accumulatation of force is only valid when GPU buffer ops are active");
224
225     GMX_ASSERT((useGpuFPmeReduction == (fPmeDeviceBuffer != nullptr)),
226                "GPU PME force reduction is only valid when a non-null GPU PME force pointer is available");
227
228     /* Skip the reduction if there was no short-range GPU work to do
229      * (either NB or both NB and bonded work). */
230     if (!pairlistIsSimple() && !haveGpuShortRangeWork(locality))
231     {
232         return;
233     }
234
235     wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
236     wallcycle_sub_start(wcycle_, ewcsNB_F_BUF_OPS);
237
238     auto fn = useGpu == BufferOpsUseGpu::True ? reduceForces<true> : reduceForces<false>;
239     fn(nbat.get(), locality, pairSearch_->gridSet(), as_rvec_array(force.data()), fPmeDeviceBuffer, pmeForcesReady, gpu_nbv, useGpuFPmeReduction, accumulateForce);
240
241     wallcycle_sub_stop(wcycle_, ewcsNB_F_BUF_OPS);
242     wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
243 }
244
245 void
246 nonbonded_verlet_t::atomdata_init_add_nbat_f_to_f_gpu()
247 {
248
249     wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
250     wallcycle_sub_start(wcycle_, ewcsNB_F_BUF_OPS);
251
252     const Nbnxm::GridSet      &gridSet = pairSearch_->gridSet();
253
254     Nbnxm::nbnxn_gpu_init_add_nbat_f_to_f(gridSet.cells().data(),
255                                           gpu_nbv,
256                                           gridSet.numRealAtomsTotal());
257
258     wallcycle_sub_stop(wcycle_, ewcsNB_F_BUF_OPS);
259     wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
260 }
261
262 real nonbonded_verlet_t::pairlistInnerRadius() const
263 {
264     return pairlistSets_->params().rlistInner;
265 }
266
267 real nonbonded_verlet_t::pairlistOuterRadius() const
268 {
269     return pairlistSets_->params().rlistOuter;
270 }
271
272 void nonbonded_verlet_t::changePairlistRadii(real rlistOuter,
273                                              real rlistInner)
274 {
275     pairlistSets_->changePairlistRadii(rlistOuter, rlistInner);
276 }
277
278 void
279 nonbonded_verlet_t::atomdata_init_copy_x_to_nbat_x_gpu()
280 {
281     Nbnxm::nbnxn_gpu_init_x_to_nbat_x(pairSearch_->gridSet(), gpu_nbv);
282 }
283
284 void nonbonded_verlet_t::insertNonlocalGpuDependency(const Nbnxm::InteractionLocality interactionLocality)
285 {
286     Nbnxm::nbnxnInsertNonlocalGpuDependency(gpu_nbv, interactionLocality);
287 }
288
289 void nonbonded_verlet_t::launch_copy_f_to_gpu(rvec *f, const Nbnxm::AtomLocality locality)
290 {
291     nbnxn_launch_copy_f_to_gpu(locality,
292                                pairSearch_->gridSet(),
293                                gpu_nbv,
294                                f);
295 }
296
297 void nonbonded_verlet_t::launch_copy_f_from_gpu(rvec *f, const Nbnxm::AtomLocality locality)
298 {
299     nbnxn_launch_copy_f_from_gpu(locality,
300                                  pairSearch_->gridSet(),
301                                  gpu_nbv,
302                                  f);
303 }
304
305 void nonbonded_verlet_t::wait_for_gpu_force_reduction(const Nbnxm::AtomLocality locality)
306 {
307     nbnxn_wait_for_gpu_force_reduction(locality, gpu_nbv);
308 }
309
310 /*! \endcond */