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