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