141b1ea9034d131224754efdc14a522a2c938e38
[alexxy/gromacs.git] / src / gromacs / nbnxm / prunekerneldispatch.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,2017,2018,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.
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 #include "gmxpre.h"
37
38 #include "gromacs/mdlib/gmx_omp_nthreads.h"
39 #include "gromacs/nbnxm/nbnxm.h"
40 #include "gromacs/timing/wallcycle.h"
41 #include "gromacs/utility/gmxassert.h"
42
43 #include "clusterdistancekerneltype.h"
44 #include "pairlistset.h"
45 #include "pairlistsets.h"
46 #include "kernels_reference/kernel_ref_prune.h"
47 #include "kernels_simd_2xmm/kernel_prune.h"
48 #include "kernels_simd_4xm/kernel_prune.h"
49
50 void PairlistSets::dispatchPruneKernel(const gmx::InteractionLocality iLocality,
51                                        const nbnxn_atomdata_t*        nbat,
52                                        const rvec*                    shift_vec)
53 {
54     pairlistSet(iLocality).dispatchPruneKernel(nbat, shift_vec);
55 }
56
57 void PairlistSet::dispatchPruneKernel(const nbnxn_atomdata_t* nbat, const rvec* shift_vec)
58 {
59     const real rlistInner = params_.rlistInner;
60
61     GMX_ASSERT(cpuLists_[0].ciOuter.size() >= cpuLists_[0].ci.size(),
62                "Here we should either have an empty ci list or ciOuter should be >= ci");
63
64     int gmx_unused nthreads = gmx_omp_nthreads_get(emntNonbonded);
65     GMX_ASSERT(nthreads == static_cast<gmx::index>(cpuLists_.size()),
66                "The number of threads should match the number of lists");
67 #pragma omp parallel for schedule(static) num_threads(nthreads)
68     for (int i = 0; i < nthreads; i++)
69     {
70         NbnxnPairlistCpu* nbl = &cpuLists_[i];
71
72         switch (getClusterDistanceKernelType(params_.pairlistType, *nbat))
73         {
74             case ClusterDistanceKernelType::CpuSimd_4xM:
75                 nbnxn_kernel_prune_4xn(nbl, nbat, shift_vec, rlistInner);
76                 break;
77             case ClusterDistanceKernelType::CpuSimd_2xMM:
78                 nbnxn_kernel_prune_2xnn(nbl, nbat, shift_vec, rlistInner);
79                 break;
80             case ClusterDistanceKernelType::CpuPlainC:
81                 nbnxn_kernel_prune_ref(nbl, nbat, shift_vec, rlistInner);
82                 break;
83             default: GMX_RELEASE_ASSERT(false, "kernel type not handled (yet)");
84         }
85     }
86 }
87
88 void nonbonded_verlet_t::dispatchPruneKernelCpu(const gmx::InteractionLocality iLocality, const rvec* shift_vec)
89 {
90     pairlistSets_->dispatchPruneKernel(iLocality, nbat.get(), shift_vec);
91 }
92
93 void nonbonded_verlet_t::dispatchPruneKernelGpu(int64_t step)
94 {
95     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
96     wallcycle_sub_start_nocount(wcycle_, ewcsLAUNCH_GPU_NONBONDED);
97
98     const bool stepIsEven = (pairlistSets().numStepsWithPairlist(step) % 2 == 0);
99
100     Nbnxm::gpu_launch_kernel_pruneonly(
101             gpu_nbv, stepIsEven ? gmx::InteractionLocality::Local : gmx::InteractionLocality::NonLocal,
102             pairlistSets().params().numRollingPruningParts);
103
104     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_NONBONDED);
105     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
106 }