3de3002b3ac3c72a6b79de09e1b59b92d8791543
[alexxy/gromacs.git] / src / gromacs / mdlib / gpuforcereduction_impl.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
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 /*! \internal \file
36  *
37  * \brief Implements backend-agnostic GPU Force Reduction functions
38  *
39  * \author Alan Gray <alang@nvidia.com>
40  *
41  * \ingroup module_mdlib
42  */
43
44 #include "gmxpre.h"
45
46 #include "gpuforcereduction_impl.h"
47
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"
53
54 namespace gmx
55 {
56
57 GpuForceReduction::Impl::Impl(const DeviceContext& deviceContext,
58                               const DeviceStream&  deviceStream,
59                               gmx_wallcycle*       wcycle) :
60     baseForce_(),
61     deviceContext_(deviceContext),
62     deviceStream_(deviceStream),
63     nbnxmForceToAdd_(),
64     rvecForceToAdd_(),
65     wcycle_(wcycle)
66 {
67 }
68
69 void GpuForceReduction::Impl::reinit(DeviceBuffer<Float3>  baseForcePtr,
70                                      const int             numAtoms,
71                                      ArrayRef<const int>   cell,
72                                      const int             atomStart,
73                                      const bool            accumulate,
74                                      GpuEventSynchronizer* completionMarker)
75 {
76     GMX_ASSERT(baseForcePtr, "Input base force for reduction has no data");
77     baseForce_        = baseForcePtr;
78     numAtoms_         = numAtoms;
79     atomStart_        = atomStart;
80     accumulate_       = accumulate;
81     completionMarker_ = completionMarker;
82     cellInfo_.cell    = cell.data();
83
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]),
89                        0,
90                        numAtoms_,
91                        deviceStream_,
92                        GpuApiCallBehavior::Async,
93                        nullptr);
94     wallcycle_stop(wcycle_, WallCycleCounter::LaunchGpu);
95
96     dependencyList_.clear();
97 };
98
99 void GpuForceReduction::Impl::registerNbnxmForce(DeviceBuffer<RVec> forcePtr)
100 {
101     GMX_ASSERT(forcePtr, "Input force for reduction has no data");
102     nbnxmForceToAdd_ = forcePtr;
103 };
104
105 void GpuForceReduction::Impl::registerRvecForce(DeviceBuffer<RVec> forcePtr)
106 {
107     GMX_ASSERT(forcePtr, "Input force for reduction has no data");
108     rvecForceToAdd_ = forcePtr;
109 };
110
111 void GpuForceReduction::Impl::addDependency(GpuEventSynchronizer* dependency)
112 {
113     dependencyList_.push_back(dependency);
114 }
115
116 void GpuForceReduction::Impl::execute()
117 {
118     wallcycle_start_nocount(wcycle_, WallCycleCounter::LaunchGpu);
119     wallcycle_sub_start(wcycle_, WallCycleSubCounter::LaunchGpuNBFBufOps);
120
121     if (numAtoms_ != 0)
122     {
123         GMX_ASSERT(nbnxmForceToAdd_, "Nbnxm force for reduction has no data");
124
125         // Enqueue wait on all dependencies passed
126         for (auto* synchronizer : dependencyList_)
127         {
128             synchronizer->enqueueWaitEvent(deviceStream_);
129         }
130
131         const bool addRvecForce = static_cast<bool>(rvecForceToAdd_); // True iff initialized
132
133         launchForceReductionKernel(numAtoms_,
134                                    atomStart_,
135                                    addRvecForce,
136                                    accumulate_,
137                                    nbnxmForceToAdd_,
138                                    rvecForceToAdd_,
139                                    baseForce_,
140                                    cellInfo_.d_cell,
141                                    deviceStream_);
142     }
143     else
144     {
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_)
150         {
151             synchronizer->consume();
152         }
153     }
154
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)
159     {
160         completionMarker_->markEvent(deviceStream_);
161     }
162
163     wallcycle_sub_stop(wcycle_, WallCycleSubCounter::LaunchGpuNBFBufOps);
164     wallcycle_stop(wcycle_, WallCycleCounter::LaunchGpu);
165 }
166
167 GpuForceReduction::GpuForceReduction(const DeviceContext& deviceContext,
168                                      const DeviceStream&  deviceStream,
169                                      gmx_wallcycle*       wcycle) :
170     impl_(new Impl(deviceContext, deviceStream, wcycle))
171 {
172 }
173
174 void GpuForceReduction::registerNbnxmForce(DeviceBuffer<RVec> forcePtr)
175 {
176     impl_->registerNbnxmForce(forcePtr);
177 }
178
179 void GpuForceReduction::registerRvecForce(DeviceBuffer<RVec> forcePtr)
180 {
181     impl_->registerRvecForce(forcePtr);
182 }
183
184 void GpuForceReduction::addDependency(GpuEventSynchronizer* dependency)
185 {
186     impl_->addDependency(dependency);
187 }
188
189 void GpuForceReduction::reinit(DeviceBuffer<RVec>    baseForcePtr,
190                                const int             numAtoms,
191                                ArrayRef<const int>   cell,
192                                const int             atomStart,
193                                const bool            accumulate,
194                                GpuEventSynchronizer* completionMarker)
195 {
196     impl_->reinit(baseForcePtr, numAtoms, cell, atomStart, accumulate, completionMarker);
197 }
198 void GpuForceReduction::execute()
199 {
200     impl_->execute();
201 }
202
203 GpuForceReduction::~GpuForceReduction() = default;
204
205 } // namespace gmx