Remove t_mdatoms from Nbnxm setAtomProperties()
[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,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.
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 "nbnxm_gpu.h"
53 #include "pairlistsets.h"
54 #include "pairsearch.h"
55
56 /*! \cond INTERNAL */
57
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)
70 {
71     nb_verlet->pairSearch_->putOnGrid(box, gridIndex, lowerCorner, upperCorner, updateGroupsCog,
72                                       atomRange, atomDensity, 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, zone, c0, c1, nullptr,
92                           { zones->cg_range[zone], zones->cg_range[zone + 1] }, -1, atomInfo, x, 0,
93                           nullptr);
94     }
95 }
96
97 bool nonbonded_verlet_t::isDynamicPruningStepCpu(int64_t step) const
98 {
99     return pairlistSets_->isDynamicPruningStepCpu(step);
100 }
101
102 bool nonbonded_verlet_t::isDynamicPruningStepGpu(int64_t step) const
103 {
104     return pairlistSets_->isDynamicPruningStepGpu(step);
105 }
106
107 gmx::ArrayRef<const int> nonbonded_verlet_t::getLocalAtomOrder() const
108 {
109     /* Return the atom order for the home cell (index 0) */
110     const Nbnxm::Grid& grid = pairSearch_->gridSet().grids()[0];
111
112     const int numIndices = grid.atomIndexEnd() - grid.firstAtomInColumn(0);
113
114     return gmx::constArrayRefFromArray(pairSearch_->gridSet().atomIndices().data(), numIndices);
115 }
116
117 void nonbonded_verlet_t::setLocalAtomOrder()
118 {
119     pairSearch_->setLocalAtomOrder();
120 }
121
122 void nonbonded_verlet_t::setAtomProperties(gmx::ArrayRef<const int>  atomTypes,
123                                            gmx::ArrayRef<const real> atomCharges,
124                                            gmx::ArrayRef<const int>  atomInfo)
125 {
126     nbnxn_atomdata_set(nbat.get(), pairSearch_->gridSet(), atomTypes, atomCharges, atomInfo);
127 }
128
129 void nonbonded_verlet_t::convertCoordinates(const gmx::AtomLocality        locality,
130                                             const bool                     fillLocal,
131                                             gmx::ArrayRef<const gmx::RVec> coordinates)
132 {
133     wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
134     wallcycle_sub_start(wcycle_, ewcsNB_X_BUF_OPS);
135
136     nbnxn_atomdata_copy_x_to_nbat_x(pairSearch_->gridSet(), locality, fillLocal,
137                                     as_rvec_array(coordinates.data()), nbat.get());
138
139     wallcycle_sub_stop(wcycle_, ewcsNB_X_BUF_OPS);
140     wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
141 }
142
143 void nonbonded_verlet_t::convertCoordinatesGpu(const gmx::AtomLocality locality,
144                                                const bool              fillLocal,
145                                                DeviceBuffer<gmx::RVec> d_x,
146                                                GpuEventSynchronizer*   xReadyOnDevice)
147 {
148     wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
149     wallcycle_sub_start(wcycle_, ewcsNB_X_BUF_OPS);
150
151     nbnxn_atomdata_x_to_nbat_x_gpu(pairSearch_->gridSet(), locality, fillLocal, gpu_nbv, d_x, xReadyOnDevice);
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 nonbonded_verlet_t::atomdata_add_nbat_f_to_f(const gmx::AtomLocality  locality,
163                                                   gmx::ArrayRef<gmx::RVec> force)
164 {
165
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))
169     {
170         return;
171     }
172
173     wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
174     wallcycle_sub_start(wcycle_, ewcsNB_F_BUF_OPS);
175
176     reduceForces(nbat.get(), locality, pairSearch_->gridSet(), as_rvec_array(force.data()));
177
178     wallcycle_sub_stop(wcycle_, ewcsNB_F_BUF_OPS);
179     wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
180 }
181
182 void nonbonded_verlet_t::atomdata_add_nbat_f_to_f_gpu(const gmx::AtomLocality locality,
183                                                       DeviceBuffer<gmx::RVec> totalForcesDevice,
184                                                       void*                   forcesPmeDevice,
185                                                       gmx::ArrayRef<GpuEventSynchronizer* const> dependencyList,
186                                                       bool useGpuFPmeReduction,
187                                                       bool accumulateForce)
188 {
189
190     GMX_ASSERT((useGpuFPmeReduction == (forcesPmeDevice != nullptr)),
191                "GPU PME force reduction is only valid when a non-null GPU PME force pointer is "
192                "available");
193
194     /* Skip the reduction if there was no short-range GPU work to do
195      * (either NB or both NB and bonded work). */
196     if (!pairlistIsSimple() && !Nbnxm::haveGpuShortRangeWork(gpu_nbv, locality))
197     {
198         return;
199     }
200
201     wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
202     wallcycle_sub_start(wcycle_, ewcsNB_F_BUF_OPS);
203
204     reduceForcesGpu(locality, totalForcesDevice, pairSearch_->gridSet(), forcesPmeDevice,
205                     dependencyList, gpu_nbv, useGpuFPmeReduction, accumulateForce);
206
207     wallcycle_sub_stop(wcycle_, ewcsNB_F_BUF_OPS);
208     wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
209 }
210
211 void nonbonded_verlet_t::atomdata_init_add_nbat_f_to_f_gpu(GpuEventSynchronizer* const localReductionDone)
212 {
213
214     wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
215     wallcycle_sub_start(wcycle_, ewcsNB_F_BUF_OPS);
216
217     const Nbnxm::GridSet& gridSet = pairSearch_->gridSet();
218
219     Nbnxm::nbnxn_gpu_init_add_nbat_f_to_f(gridSet.cells().data(), gpu_nbv,
220                                           gridSet.numRealAtomsTotal(), localReductionDone);
221
222     wallcycle_sub_stop(wcycle_, ewcsNB_F_BUF_OPS);
223     wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
224 }
225
226 real nonbonded_verlet_t::pairlistInnerRadius() const
227 {
228     return pairlistSets_->params().rlistInner;
229 }
230
231 real nonbonded_verlet_t::pairlistOuterRadius() const
232 {
233     return pairlistSets_->params().rlistOuter;
234 }
235
236 void nonbonded_verlet_t::changePairlistRadii(real rlistOuter, real rlistInner)
237 {
238     pairlistSets_->changePairlistRadii(rlistOuter, rlistInner);
239 }
240
241 void nonbonded_verlet_t::setupGpuShortRangeWork(const gmx::GpuBonded*          gpuBonded,
242                                                 const gmx::InteractionLocality iLocality)
243 {
244     if (useGpu() && !emulateGpu())
245     {
246         Nbnxm::setupGpuShortRangeWork(gpu_nbv, gpuBonded, iLocality);
247     }
248 }
249
250 void nonbonded_verlet_t::atomdata_init_copy_x_to_nbat_x_gpu()
251 {
252     Nbnxm::nbnxn_gpu_init_x_to_nbat_x(pairSearch_->gridSet(), gpu_nbv);
253 }
254
255 void nonbonded_verlet_t::insertNonlocalGpuDependency(const gmx::InteractionLocality interactionLocality)
256 {
257     Nbnxm::nbnxnInsertNonlocalGpuDependency(gpu_nbv, interactionLocality);
258 }
259
260 /*! \endcond */