2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2020,2021, 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.
37 * \brief Implements backend-agnostic GPU Force Reduction functions
39 * \author Alan Gray <alang@nvidia.com>
41 * \ingroup module_mdlib
46 #include "gpuforcereduction_impl.h"
48 #include "gromacs/gpu_utils/device_stream.h"
49 #include "gromacs/gpu_utils/devicebuffer.h"
50 #include "gromacs/gpu_utils/gpueventsynchronizer.h"
51 #include "gromacs/mdlib/gpuforcereduction_impl_internal.h"
52 #include "gromacs/utility/gmxassert.h"
57 GpuForceReduction::Impl::Impl(const DeviceContext& deviceContext,
58 const DeviceStream& deviceStream,
59 gmx_wallcycle* wcycle) :
61 deviceContext_(deviceContext),
62 deviceStream_(deviceStream),
69 void GpuForceReduction::Impl::reinit(DeviceBuffer<Float3> baseForcePtr,
71 ArrayRef<const int> cell,
73 const bool accumulate,
74 GpuEventSynchronizer* completionMarker)
76 GMX_ASSERT(baseForcePtr, "Input base force for reduction has no data");
77 baseForce_ = baseForcePtr;
79 atomStart_ = atomStart;
80 accumulate_ = accumulate;
81 completionMarker_ = completionMarker;
82 cellInfo_.cell = cell.data();
84 wallcycle_start_nocount(wcycle_, WallCycleCounter::LaunchGpu);
85 reallocateDeviceBuffer(
86 &cellInfo_.d_cell, numAtoms_, &cellInfo_.cellSize, &cellInfo_.cellSizeAlloc, deviceContext_);
87 copyToDeviceBuffer(&cellInfo_.d_cell,
88 &(cellInfo_.cell[atomStart]),
92 GpuApiCallBehavior::Async,
94 wallcycle_stop(wcycle_, WallCycleCounter::LaunchGpu);
96 dependencyList_.clear();
99 void GpuForceReduction::Impl::registerNbnxmForce(DeviceBuffer<RVec> forcePtr)
101 GMX_ASSERT(forcePtr, "Input force for reduction has no data");
102 nbnxmForceToAdd_ = forcePtr;
105 void GpuForceReduction::Impl::registerRvecForce(DeviceBuffer<RVec> forcePtr)
107 GMX_ASSERT(forcePtr, "Input force for reduction has no data");
108 rvecForceToAdd_ = forcePtr;
111 void GpuForceReduction::Impl::addDependency(GpuEventSynchronizer* dependency)
113 dependencyList_.push_back(dependency);
116 void GpuForceReduction::Impl::execute()
118 wallcycle_start_nocount(wcycle_, WallCycleCounter::LaunchGpu);
119 wallcycle_sub_start(wcycle_, WallCycleSubCounter::LaunchGpuNBFBufOps);
123 GMX_ASSERT(nbnxmForceToAdd_, "Nbnxm force for reduction has no data");
125 // Enqueue wait on all dependencies passed
126 for (auto* synchronizer : dependencyList_)
128 synchronizer->enqueueWaitEvent(deviceStream_);
131 const bool addRvecForce = static_cast<bool>(rvecForceToAdd_); // True iff initialized
133 launchForceReductionKernel(numAtoms_,
145 /* In case we have nothing to do, but still have dependencies, we need
146 * to consume them and mark our own event.
147 * Happens sometimes in MdrunVsitesTest.
148 * Issue #3988, #4227. */
149 for (auto* synchronizer : dependencyList_)
151 synchronizer->consume();
155 /* Mark that kernel has been launched.
156 * Even if we have no work to do and have not launched the kernel, we still mark the event
157 * in order to ensure proper marking/consumption balance, see Issue #3988, #4227. */
158 if (completionMarker_ != nullptr)
160 completionMarker_->markEvent(deviceStream_);
163 wallcycle_sub_stop(wcycle_, WallCycleSubCounter::LaunchGpuNBFBufOps);
164 wallcycle_stop(wcycle_, WallCycleCounter::LaunchGpu);
167 GpuForceReduction::GpuForceReduction(const DeviceContext& deviceContext,
168 const DeviceStream& deviceStream,
169 gmx_wallcycle* wcycle) :
170 impl_(new Impl(deviceContext, deviceStream, wcycle))
174 void GpuForceReduction::registerNbnxmForce(DeviceBuffer<RVec> forcePtr)
176 impl_->registerNbnxmForce(forcePtr);
179 void GpuForceReduction::registerRvecForce(DeviceBuffer<RVec> forcePtr)
181 impl_->registerRvecForce(forcePtr);
184 void GpuForceReduction::addDependency(GpuEventSynchronizer* dependency)
186 impl_->addDependency(dependency);
189 void GpuForceReduction::reinit(DeviceBuffer<RVec> baseForcePtr,
191 ArrayRef<const int> cell,
193 const bool accumulate,
194 GpuEventSynchronizer* completionMarker)
196 impl_->reinit(baseForcePtr, numAtoms, cell, atomStart, accumulate, completionMarker);
198 void GpuForceReduction::execute()
203 GpuForceReduction::~GpuForceReduction() = default;