00a15a0c38412e0f629c84043256bd8cc0699932
[alexxy/gromacs.git] / src / gromacs / nbnxm / nbnxm_gpu_buffer_ops.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
36 /*! \internal \file
37  *  \brief
38  *  Common code for GPU buffer operations, namely the coordinate layout conversion
39  *
40  *  \ingroup module_nbnxm
41  */
42 #include "gmxpre.h"
43
44 #include "config.h"
45
46 #include "gromacs/gpu_utils/device_stream.h"
47 #if GMX_GPU_CUDA
48 #    include "gromacs/gpu_utils/gpueventsynchronizer.cuh"
49 #endif
50 #include "gromacs/mdtypes/locality.h"
51 #include "gromacs/nbnxm/gridset.h"
52 #include "gromacs/nbnxm/nbnxm_gpu.h"
53 #include "gromacs/nbnxm/nbnxm_gpu_buffer_ops_internal.h"
54 #if GMX_GPU_CUDA
55 #    include "gromacs/nbnxm/cuda/nbnxm_cuda_types.h"
56 #endif
57 #include "gromacs/utility/exceptions.h"
58
59 namespace Nbnxm
60 {
61
62 void nbnxn_gpu_x_to_nbat_x(const Nbnxm::Grid&      grid,
63                            NbnxmGpu*               nb,
64                            DeviceBuffer<gmx::RVec> d_x,
65                            GpuEventSynchronizer*   xReadyOnDevice,
66                            const gmx::AtomLocality locality,
67                            int                     gridId,
68                            int                     numColumnsMax,
69                            bool                    mustInsertNonLocalDependency)
70 {
71     GMX_RELEASE_ASSERT(GMX_GPU_CUDA, "nbnxn_gpu_x_to_nbat_x only supported with CUDA");
72     GMX_ASSERT(nb, "Need a valid nbnxn_gpu object");
73     gmx::InteractionLocality interactionLoc = gmx::atomToInteractionLocality(locality);
74
75     const DeviceStream& deviceStream = *nb->deviceStreams[interactionLoc];
76
77     const int numAtoms = grid.srcAtomEnd() - grid.srcAtomBegin();
78
79     // Only insert wait on the first iteration of the loop.
80     if (xReadyOnDevice != nullptr)
81     {
82         xReadyOnDevice->enqueueWaitEvent(deviceStream);
83     }
84
85     // avoid empty kernel launch, skip to inserting stream dependency
86     if (numAtoms != 0)
87     {
88         GMX_ASSERT(d_x, "Need a valid device pointer");
89         launchNbnxmKernelTransformXToXq(grid, nb, d_x, deviceStream, numColumnsMax, gridId);
90     }
91
92     if (mustInsertNonLocalDependency)
93     {
94         Nbnxm::nbnxnInsertNonlocalGpuDependency(nb, interactionLoc);
95     }
96 }
97
98 DeviceBuffer<Float3> getGpuForces(NbnxmGpu* nb)
99 {
100     return nb->atdat->f;
101 }
102
103 } // namespace Nbnxm