Fixes for F buffer ops change
[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::setCoordinates(const Nbnxm::AtomLocality       locality,
136                                         const bool                      fillLocal,
137                                         gmx::ArrayRef<const gmx::RVec>  x,
138                                         BufferOpsUseGpu                 useGpu,
139                                         void                           *xPmeDevicePtr,
140                                         gmx_wallcycle                  *wcycle)
141 {
142     wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
143     wallcycle_sub_start(wcycle, ewcsNB_X_BUF_OPS);
144
145     auto fnPtr = (useGpu == BufferOpsUseGpu::True) ?
146         nbnxn_atomdata_copy_x_to_nbat_x<true> :
147         nbnxn_atomdata_copy_x_to_nbat_x<false>;
148
149     fnPtr(pairSearch_->gridSet(), locality, fillLocal,
150           as_rvec_array(x.data()),
151           nbat.get(), gpu_nbv, xPmeDevicePtr);
152
153     wallcycle_sub_stop(wcycle, ewcsNB_X_BUF_OPS);
154     wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
155 }
156
157 gmx::ArrayRef<const int> nonbonded_verlet_t::getGridIndices() const
158 {
159     return pairSearch_->gridSet().cells();
160 }
161
162 void
163 nonbonded_verlet_t::atomdata_add_nbat_f_to_f(const Nbnxm::AtomLocality           locality,
164                                              rvec                               *f,
165                                              BufferOpsUseGpu                     useGpu,
166                                              GpuBufferOpsAccumulateForce         accumulateForce,
167                                              gmx_wallcycle                      *wcycle)
168 {
169
170     GMX_ASSERT(!((useGpu == BufferOpsUseGpu::False) &&
171                  (accumulateForce == GpuBufferOpsAccumulateForce::True)),
172                "Accumulatation of force is only valid when GPU buffer ops are active");
173
174     /* Skip the reduction if there was no short-range GPU work to do
175      * (either NB or both NB and bonded work). */
176     if (!pairlistIsSimple() && !haveGpuShortRangeWork(locality))
177     {
178         return;
179     }
180
181     wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
182     wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
183
184     auto fn = useGpu == BufferOpsUseGpu::True ? reduceForces<true> : reduceForces<false>;
185     fn(nbat.get(), locality, pairSearch_->gridSet(), f, gpu_nbv, accumulateForce);
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_init_add_nbat_f_to_f_gpu(gmx_wallcycle *wcycle)
193 {
194
195     wallcycle_start(wcycle, ewcNB_XF_BUF_OPS);
196     wallcycle_sub_start(wcycle, ewcsNB_F_BUF_OPS);
197
198     const Nbnxm::GridSet      &gridSet = pairSearch_->gridSet();
199
200     Nbnxm::nbnxn_gpu_init_add_nbat_f_to_f(gridSet.cells().data(),
201                                           gpu_nbv,
202                                           gridSet.numRealAtomsTotal());
203
204     wallcycle_sub_stop(wcycle, ewcsNB_F_BUF_OPS);
205     wallcycle_stop(wcycle, ewcNB_XF_BUF_OPS);
206 }
207
208 real nonbonded_verlet_t::pairlistInnerRadius() const
209 {
210     return pairlistSets_->params().rlistInner;
211 }
212
213 real nonbonded_verlet_t::pairlistOuterRadius() const
214 {
215     return pairlistSets_->params().rlistOuter;
216 }
217
218 void nonbonded_verlet_t::changePairlistRadii(real rlistOuter,
219                                              real rlistInner)
220 {
221     pairlistSets_->changePairlistRadii(rlistOuter, rlistInner);
222 }
223
224 void
225 nonbonded_verlet_t::atomdata_init_copy_x_to_nbat_x_gpu()
226 {
227     Nbnxm::nbnxn_gpu_init_x_to_nbat_x(pairSearch_->gridSet(), gpu_nbv);
228 }
229
230 void nonbonded_verlet_t::insertNonlocalGpuDependency(const Nbnxm::InteractionLocality interactionLocality)
231 {
232     Nbnxm::nbnxnInsertNonlocalGpuDependency(gpu_nbv, interactionLocality);
233 }
234
235 void nonbonded_verlet_t::launch_copy_f_to_gpu(rvec *f, const Nbnxm::AtomLocality locality)
236 {
237     nbnxn_launch_copy_f_to_gpu(locality,
238                                pairSearch_->gridSet(),
239                                gpu_nbv,
240                                f);
241 }
242
243 void nonbonded_verlet_t::launch_copy_f_from_gpu(rvec *f, const Nbnxm::AtomLocality locality)
244 {
245     nbnxn_launch_copy_f_from_gpu(locality,
246                                  pairSearch_->gridSet(),
247                                  gpu_nbv,
248                                  f);
249 }
250
251 void nonbonded_verlet_t::wait_for_gpu_force_reduction(const Nbnxm::AtomLocality locality)
252 {
253     nbnxn_wait_for_gpu_force_reduction(locality, gpu_nbv);
254 }
255
256 /*! \endcond */