Introduce SimulatorBuilder
[gromacs.git] / src / gromacs / nbnxm / nbnxm.cpp
blob77119044cbcd08cc30362feb737b8c6e5ce4d0ea
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
36 /*! \internal \file
37 * \brief
38 * Implements the Nbnxm class
40 * \author Berk Hess <hess@kth.se>
41 * \ingroup module_nbnxm
44 #include "gmxpre.h"
46 #include "nbnxm.h"
48 #include "gromacs/domdec/domdec_struct.h"
49 #include "gromacs/nbnxm/atomdata.h"
50 #include "gromacs/timing/wallcycle.h"
52 #include "pairlistsets.h"
53 #include "pairsearch.h"
55 /*! \cond INTERNAL */
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)
71 nb_verlet->pairSearch_->putOnGrid(box, ddZone, lowerCorner, upperCorner,
72 updateGroupsCog, atomStart, atomEnd, atomDensity,
73 atomInfo, x, numAtomsMoved, move,
74 nb_verlet->nbat.get());
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)
83 for (int zone = 1; zone < zones->n; zone++)
85 rvec c0, c1;
86 for (int d = 0; d < DIM; d++)
88 c0[d] = zones->size[zone].bb_x0[d];
89 c1[d] = zones->size[zone].bb_x1[d];
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,
100 0, nullptr);
104 bool nonbonded_verlet_t::isDynamicPruningStepCpu(int64_t step) const
106 return pairlistSets_->isDynamicPruningStepCpu(step);
109 bool nonbonded_verlet_t::isDynamicPruningStepGpu(int64_t step) const
111 return pairlistSets_->isDynamicPruningStepGpu(step);
114 gmx::ArrayRef<const int> nonbonded_verlet_t::getLocalAtomOrder() const
116 /* Return the atom order for the home cell (index 0) */
117 const Nbnxm::Grid &grid = pairSearch_->gridSet().grids()[0];
119 const int numIndices = grid.atomIndexEnd() - grid.firstAtomInColumn(0);
121 return gmx::constArrayRefFromArray(pairSearch_->gridSet().atomIndices().data(), numIndices);
124 void nonbonded_verlet_t::setLocalAtomOrder()
126 pairSearch_->setLocalAtomOrder();
129 void nonbonded_verlet_t::setAtomProperties(const t_mdatoms &mdatoms,
130 gmx::ArrayRef<const int> atomInfo)
132 nbnxn_atomdata_set(nbat.get(), pairSearch_->gridSet(), &mdatoms, atomInfo.data());
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)
141 wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
142 wallcycle_sub_start(wcycle_, ewcsNB_X_BUF_OPS);
144 auto fnPtr = (useGpu == BufferOpsUseGpu::True) ?
145 nbnxn_atomdata_copy_x_to_nbat_x<true> :
146 nbnxn_atomdata_copy_x_to_nbat_x<false>;
148 fnPtr(pairSearch_->gridSet(), locality, fillLocal,
149 as_rvec_array(x.data()),
150 nbat.get(), gpu_nbv, xPmeDevicePtr);
152 wallcycle_sub_stop(wcycle_, ewcsNB_X_BUF_OPS);
153 wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
156 gmx::ArrayRef<const int> nonbonded_verlet_t::getGridIndices() const
158 return pairSearch_->gridSet().cells();
161 void
162 nonbonded_verlet_t::atomdata_add_nbat_f_to_f(const Nbnxm::AtomLocality locality,
163 rvec *f,
164 BufferOpsUseGpu useGpu,
165 GpuBufferOpsAccumulateForce accumulateForce)
168 GMX_ASSERT(!((useGpu == BufferOpsUseGpu::False) &&
169 (accumulateForce == GpuBufferOpsAccumulateForce::True)),
170 "Accumulatation of force is only valid when GPU buffer ops are active");
172 /* Skip the reduction if there was no short-range GPU work to do
173 * (either NB or both NB and bonded work). */
174 if (!pairlistIsSimple() && !haveGpuShortRangeWork(locality))
176 return;
179 wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
180 wallcycle_sub_start(wcycle_, ewcsNB_F_BUF_OPS);
182 auto fn = useGpu == BufferOpsUseGpu::True ? reduceForces<true> : reduceForces<false>;
183 fn(nbat.get(), locality, pairSearch_->gridSet(), f, gpu_nbv, accumulateForce);
185 wallcycle_sub_stop(wcycle_, ewcsNB_F_BUF_OPS);
186 wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
189 void
190 nonbonded_verlet_t::atomdata_init_add_nbat_f_to_f_gpu()
193 wallcycle_start(wcycle_, ewcNB_XF_BUF_OPS);
194 wallcycle_sub_start(wcycle_, ewcsNB_F_BUF_OPS);
196 const Nbnxm::GridSet &gridSet = pairSearch_->gridSet();
198 Nbnxm::nbnxn_gpu_init_add_nbat_f_to_f(gridSet.cells().data(),
199 gpu_nbv,
200 gridSet.numRealAtomsTotal());
202 wallcycle_sub_stop(wcycle_, ewcsNB_F_BUF_OPS);
203 wallcycle_stop(wcycle_, ewcNB_XF_BUF_OPS);
206 real nonbonded_verlet_t::pairlistInnerRadius() const
208 return pairlistSets_->params().rlistInner;
211 real nonbonded_verlet_t::pairlistOuterRadius() const
213 return pairlistSets_->params().rlistOuter;
216 void nonbonded_verlet_t::changePairlistRadii(real rlistOuter,
217 real rlistInner)
219 pairlistSets_->changePairlistRadii(rlistOuter, rlistInner);
222 void
223 nonbonded_verlet_t::atomdata_init_copy_x_to_nbat_x_gpu()
225 Nbnxm::nbnxn_gpu_init_x_to_nbat_x(pairSearch_->gridSet(), gpu_nbv);
228 void nonbonded_verlet_t::insertNonlocalGpuDependency(const Nbnxm::InteractionLocality interactionLocality)
230 Nbnxm::nbnxnInsertNonlocalGpuDependency(gpu_nbv, interactionLocality);
233 void nonbonded_verlet_t::launch_copy_f_to_gpu(rvec *f, const Nbnxm::AtomLocality locality)
235 nbnxn_launch_copy_f_to_gpu(locality,
236 pairSearch_->gridSet(),
237 gpu_nbv,
241 void nonbonded_verlet_t::launch_copy_f_from_gpu(rvec *f, const Nbnxm::AtomLocality locality)
243 nbnxn_launch_copy_f_from_gpu(locality,
244 pairSearch_->gridSet(),
245 gpu_nbv,
249 void nonbonded_verlet_t::wait_for_gpu_force_reduction(const Nbnxm::AtomLocality locality)
251 nbnxn_wait_for_gpu_force_reduction(locality, gpu_nbv);
254 /*! \endcond */