2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2017, 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 * \brief Declare interface for GPU execution for NBNXN module
38 * \author Szilard Pall <pall.szilard@gmail.com>
39 * \author Mark Abraham <mark.j.abraham@gmail.com>
40 * \ingroup module_mdlib
43 #ifndef GMX_MDLIB_NBNXN_GPU_H
44 #define GMX_MDLIB_NBNXN_GPU_H
46 #include "gromacs/gpu_utils/gpu_macros.h"
47 #include "gromacs/math/vectypes.h"
48 #include "gromacs/mdlib/nbnxn_gpu_types.h"
49 #include "gromacs/utility/basedefinitions.h"
50 #include "gromacs/utility/real.h"
56 struct nbnxn_atomdata_t;
59 * Launch asynchronously the nonbonded force calculations.
61 * This consists of the following (async) steps launched:
63 * - upload shift vector;
65 * The local and non-local interaction calculations are launched in two
69 void nbnxn_gpu_launch_kernel(gmx_nbnxn_gpu_t gmx_unused *nb,
70 const struct nbnxn_atomdata_t gmx_unused *nbdata,
72 int gmx_unused iloc) GPU_FUNC_TERM
75 * Launch asynchronously the nonbonded prune-only kernel.
77 * The local and non-local list pruning are launched in their separate streams.
79 * Notes for future scheduling tuning:
80 * Currently we schedule the dynamic pruning between two MD steps *after* both local and
81 * nonlocal force D2H transfers completed. We could launch already after the cpyback
82 * is launched, but we want to avoid prune kernels (especially in the non-local
83 * high prio-stream) competing with nonbonded work.
85 * However, this is not ideal as this schedule does not expose the available
86 * concurrency. The dynamic pruning kernel:
87 * - should be allowed to overlap with any task other than force compute, including
88 * transfers (F D2H and the next step's x H2D as well as force clearing).
89 * - we'd prefer to avoid competition with non-bonded force kernels belonging
90 * to the same rank and ideally other ranks too.
92 * In the most general case, the former would require scheduling pruning in a separate
93 * stream and adding additional event sync points to ensure that force kernels read
94 * consistent pair list data. This would lead to some overhead (due to extra
95 * cudaStreamWaitEvent calls, 3-5 us/call) which we might be able to live with.
96 * The gains from additional overlap might not be significant as long as
97 * update+constraints anyway takes longer than pruning, but there will still
98 * be use-cases where more overlap may help (e.g. multiple ranks per GPU,
99 * no/hbonds only constraints).
100 * The above second point is harder to address given that multiple ranks will often
101 * share a GPU. Ranks that complete their nonbondeds sooner can schedule pruning earlier
102 * and without a third priority level it is difficult to avoid some interference of
103 * prune kernels with force tasks (in particular preemption of low-prio local force task).
105 * \param [inout] nb GPU nonbonded data.
106 * \param [in] iloc Interaction locality flag.
107 * \param [in] numParts Number of parts the pair list is split into in the rolling kernel.
110 void nbnxn_gpu_launch_kernel_pruneonly(gmx_nbnxn_gpu_t gmx_unused *nb,
112 int gmx_unused numParts) GPU_FUNC_TERM
115 * Launch asynchronously the download of nonbonded forces from the GPU
116 * (and energies/shift forces if required).
119 void nbnxn_gpu_launch_cpyback(gmx_nbnxn_gpu_t gmx_unused *nb,
120 const struct nbnxn_atomdata_t gmx_unused *nbatom,
121 int gmx_unused flags,
122 int gmx_unused aloc) GPU_FUNC_TERM
125 * Wait for the asynchronously launched nonbonded calculations and data
126 * transfers to finish.
129 void nbnxn_gpu_wait_for_gpu(gmx_nbnxn_gpu_t gmx_unused *nb,
130 int gmx_unused flags,
132 real gmx_unused *e_lj,
133 real gmx_unused *e_el,
134 rvec gmx_unused *fshift) GPU_FUNC_TERM
136 /*! \brief Selects the Ewald kernel type, analytical or tabulated, single or twin cut-off. */
138 int nbnxn_gpu_pick_ewald_kernel_type(bool gmx_unused bTwinCut) GPU_FUNC_TERM_WITH_RETURN(-1)