2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
5 * Copyright (c) 2017,2018,2019,2020,2021, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
39 #include "gromacs/gmxlib/nrnb.h"
40 #include "gromacs/math/vectypes.h"
41 #include "gromacs/mdlib/enerdata_utils.h"
42 #include "gromacs/mdlib/force.h"
43 #include "gromacs/mdlib/gmx_omp_nthreads.h"
44 #include "gromacs/mdtypes/enerdata.h"
45 #include "gromacs/mdtypes/forceoutput.h"
46 #include "gromacs/mdtypes/forcerec.h"
47 #include "gromacs/mdtypes/inputrec.h"
48 #include "gromacs/mdtypes/interaction_const.h"
49 #include "gromacs/mdtypes/md_enums.h"
50 #include "gromacs/mdtypes/mdatom.h"
51 #include "gromacs/mdtypes/nblist.h"
52 #include "gromacs/mdtypes/simulation_workload.h"
53 #include "gromacs/nbnxm/gpu_data_mgmt.h"
54 #include "gromacs/nbnxm/nbnxm.h"
55 #include "gromacs/simd/simd.h"
56 #include "gromacs/timing/wallcycle.h"
57 #include "gromacs/utility/enumerationhelpers.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/gmxassert.h"
60 #include "gromacs/utility/real.h"
62 #include "kernel_common.h"
63 #include "nbnxm_gpu.h"
64 #include "nbnxm_simd.h"
65 #include "pairlistset.h"
66 #include "pairlistsets.h"
67 #include "kernels_reference/kernel_gpu_ref.h"
68 #define INCLUDE_KERNELFUNCTION_TABLES
69 #include "kernels_reference/kernel_ref.h"
70 #ifdef GMX_NBNXN_SIMD_2XNN
71 # include "kernels_simd_2xmm/kernels.h"
73 #ifdef GMX_NBNXN_SIMD_4XN
74 # include "kernels_simd_4xm/kernels.h"
76 #undef INCLUDE_FUNCTION_TABLES
78 /*! \brief Clears the energy group output buffers
80 * \param[in,out] out nbnxn kernel output struct
82 static void clearGroupEnergies(nbnxn_atomdata_output_t* out)
84 std::fill(out->Vvdw.begin(), out->Vvdw.end(), 0.0_real);
85 std::fill(out->Vc.begin(), out->Vc.end(), 0.0_real);
86 std::fill(out->VSvdw.begin(), out->VSvdw.end(), 0.0_real);
87 std::fill(out->VSc.begin(), out->VSc.end(), 0.0_real);
90 /*! \brief Reduce the group-pair energy buffers produced by a SIMD kernel
91 * to single terms in the output buffers.
93 * The SIMD kernels produce a large number of energy buffer in SIMD registers
94 * to avoid scattered reads and writes.
96 * \tparam unrollj The unroll size for j-particles in the SIMD kernel
97 * \param[in] numGroups The number of energy groups
98 * \param[in] numGroups_2log Log2 of numGroups, rounded up
99 * \param[in,out] out Struct with energy buffers
101 template<int unrollj>
102 static void reduceGroupEnergySimdBuffers(int numGroups, int numGroups_2log, nbnxn_atomdata_output_t* out)
104 const int unrollj_half = unrollj / 2;
105 /* Energies are stored in SIMD registers with size 2^numGroups_2log */
106 const int numGroupsStorage = (1 << numGroups_2log);
108 const real* gmx_restrict vVdwSimd = out->VSvdw.data();
109 const real* gmx_restrict vCoulombSimd = out->VSc.data();
110 real* gmx_restrict vVdw = out->Vvdw.data();
111 real* gmx_restrict vCoulomb = out->Vc.data();
113 /* The size of the SIMD energy group buffer array is:
114 * numGroups*numGroups*numGroupsStorage*unrollj_half*simd_width
116 for (int i = 0; i < numGroups; i++)
118 for (int j1 = 0; j1 < numGroups; j1++)
120 for (int j0 = 0; j0 < numGroups; j0++)
122 int c = ((i * numGroups + j1) * numGroupsStorage + j0) * unrollj_half * unrollj;
123 for (int s = 0; s < unrollj_half; s++)
125 vVdw[i * numGroups + j0] += vVdwSimd[c + 0];
126 vVdw[i * numGroups + j1] += vVdwSimd[c + 1];
127 vCoulomb[i * numGroups + j0] += vCoulombSimd[c + 0];
128 vCoulomb[i * numGroups + j1] += vCoulombSimd[c + 1];
136 CoulombKernelType getCoulombKernelType(const Nbnxm::EwaldExclusionType ewaldExclusionType,
137 const CoulombInteractionType coulombInteractionType,
138 const bool haveEqualCoulombVwdRadii)
141 if (EEL_RF(coulombInteractionType) || coulombInteractionType == CoulombInteractionType::Cut)
143 return CoulombKernelType::ReactionField;
147 if (ewaldExclusionType == Nbnxm::EwaldExclusionType::Table)
149 if (haveEqualCoulombVwdRadii)
151 return CoulombKernelType::Table;
155 return CoulombKernelType::TableTwin;
160 if (haveEqualCoulombVwdRadii)
162 return CoulombKernelType::Ewald;
166 return CoulombKernelType::EwaldTwin;
172 int getVdwKernelType(const Nbnxm::KernelType kernelType,
173 const LJCombinationRule ljCombinationRule,
174 const VanDerWaalsType vanDerWaalsType,
175 const InteractionModifiers interactionModifiers,
176 const LongRangeVdW longRangeVdW)
178 if (vanDerWaalsType == VanDerWaalsType::Cut)
180 switch (interactionModifiers)
182 case InteractionModifiers::None:
183 case InteractionModifiers::PotShift:
184 switch (ljCombinationRule)
186 case LJCombinationRule::Geometric: return vdwktLJCUT_COMBGEOM;
187 case LJCombinationRule::LorentzBerthelot: return vdwktLJCUT_COMBLB;
188 case LJCombinationRule::None: return vdwktLJCUT_COMBNONE;
189 default: GMX_THROW(gmx::InvalidInputError("Unknown combination rule"));
191 case InteractionModifiers::ForceSwitch: return vdwktLJFORCESWITCH;
192 case InteractionModifiers::PotSwitch: return vdwktLJPOTSWITCH;
194 std::string errorMsg =
195 gmx::formatString("Unsupported VdW interaction modifier %s (%d)",
196 enumValueToString(interactionModifiers),
197 static_cast<int>(interactionModifiers));
198 GMX_THROW(gmx::InvalidInputError(errorMsg));
201 else if (vanDerWaalsType == VanDerWaalsType::Pme)
203 if (longRangeVdW == LongRangeVdW::Geom)
205 return vdwktLJEWALDCOMBGEOM;
209 /* At setup we (should have) selected the C reference kernel */
210 GMX_RELEASE_ASSERT(kernelType == Nbnxm::KernelType::Cpu4x4_PlainC,
211 "Only the C reference nbnxn SIMD kernel supports LJ-PME with LB "
212 "combination rules");
213 return vdwktLJEWALDCOMBLB;
218 std::string errorMsg = gmx::formatString("Unsupported VdW interaction type %s (%d)",
219 enumValueToString(vanDerWaalsType),
220 static_cast<int>(vanDerWaalsType));
221 GMX_THROW(gmx::InvalidInputError(errorMsg));
225 /*! \brief Dispatches the non-bonded N versus M atom cluster CPU kernels.
227 * OpenMP parallelization is performed within this function.
228 * Energy reduction, but not force and shift force reduction, is performed
229 * within this function.
231 * \param[in] pairlistSet Pairlists with local or non-local interactions to compute
232 * \param[in] kernelSetup The non-bonded kernel setup
233 * \param[in,out] nbat The atomdata for the interactions
234 * \param[in] ic Non-bonded interaction constants
235 * \param[in] shiftVectors The PBC shift vectors
236 * \param[in] stepWork Flags that tell what to compute
237 * \param[in] clearF Enum that tells if to clear the force output buffer
238 * \param[out] vCoulomb Output buffer for Coulomb energies
239 * \param[out] vVdw Output buffer for Van der Waals energies
240 * \param[in] wcycle Pointer to cycle counting data structure.
242 static void nbnxn_kernel_cpu(const PairlistSet& pairlistSet,
243 const Nbnxm::KernelSetup& kernelSetup,
244 nbnxn_atomdata_t* nbat,
245 const interaction_const_t& ic,
246 gmx::ArrayRef<const gmx::RVec> shiftVectors,
247 const gmx::StepWorkload& stepWork,
251 gmx_wallcycle* wcycle)
254 const nbnxn_atomdata_t::Params& nbatParams = nbat->params();
256 const int coulkt = static_cast<int>(getCoulombKernelType(
257 kernelSetup.ewaldExclusionType, ic.eeltype, (ic.rcoulomb == ic.rvdw)));
258 const int vdwkt = getVdwKernelType(
259 kernelSetup.kernelType, nbatParams.ljCombinationRule, ic.vdwtype, ic.vdw_modifier, ic.ljpme_comb_rule);
261 gmx::ArrayRef<const NbnxnPairlistCpu> pairlists = pairlistSet.cpuLists();
263 const auto* shiftVecPointer = as_rvec_array(shiftVectors.data());
265 int gmx_unused nthreads = gmx_omp_nthreads_get(ModuleMultiThread::Nonbonded);
266 wallcycle_sub_start(wcycle, WallCycleSubCounter::NonbondedClear);
267 #pragma omp parallel for schedule(static) num_threads(nthreads)
268 for (gmx::index nb = 0; nb < pairlists.ssize(); nb++)
270 // Presently, the kernels do not call C++ code that can throw,
271 // so no need for a try/catch pair in this OpenMP region.
272 nbnxn_atomdata_output_t* out = &nbat->out[nb];
274 if (clearF == enbvClearFYes)
276 clearForceBuffer(nbat, nb);
278 clear_fshift(out->fshift.data());
283 wallcycle_sub_stop(wcycle, WallCycleSubCounter::NonbondedClear);
284 wallcycle_sub_start(wcycle, WallCycleSubCounter::NonbondedKernel);
287 // TODO: Change to reference
288 const NbnxnPairlistCpu* pairlist = &pairlists[nb];
290 if (!stepWork.computeEnergy)
292 /* Don't calculate energies */
293 switch (kernelSetup.kernelType)
295 case Nbnxm::KernelType::Cpu4x4_PlainC:
296 nbnxn_kernel_noener_ref[coulkt][vdwkt](pairlist, nbat, &ic, shiftVecPointer, out);
298 #ifdef GMX_NBNXN_SIMD_2XNN
299 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
300 nbnxm_kernel_noener_simd_2xmm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVecPointer, out);
303 #ifdef GMX_NBNXN_SIMD_4XN
304 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
305 nbnxm_kernel_noener_simd_4xm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVecPointer, out);
308 default: GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
311 else if (out->Vvdw.size() == 1)
313 /* A single energy group (pair) */
317 switch (kernelSetup.kernelType)
319 case Nbnxm::KernelType::Cpu4x4_PlainC:
320 nbnxn_kernel_ener_ref[coulkt][vdwkt](pairlist, nbat, &ic, shiftVecPointer, out);
322 #ifdef GMX_NBNXN_SIMD_2XNN
323 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
324 nbnxm_kernel_ener_simd_2xmm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVecPointer, out);
327 #ifdef GMX_NBNXN_SIMD_4XN
328 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
329 nbnxm_kernel_ener_simd_4xm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVecPointer, out);
332 default: GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
337 /* Calculate energy group contributions */
338 clearGroupEnergies(out);
342 switch (kernelSetup.kernelType)
344 case Nbnxm::KernelType::Cpu4x4_PlainC:
345 unrollj = c_nbnxnCpuIClusterSize;
346 nbnxn_kernel_energrp_ref[coulkt][vdwkt](pairlist, nbat, &ic, shiftVecPointer, out);
348 #ifdef GMX_NBNXN_SIMD_2XNN
349 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
350 unrollj = GMX_SIMD_REAL_WIDTH / 2;
351 nbnxm_kernel_energrp_simd_2xmm[coulkt][vdwkt](
352 pairlist, nbat, &ic, shiftVecPointer, out);
355 #ifdef GMX_NBNXN_SIMD_4XN
356 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
357 unrollj = GMX_SIMD_REAL_WIDTH;
358 nbnxm_kernel_energrp_simd_4xm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVecPointer, out);
361 default: GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
364 if (kernelSetup.kernelType != Nbnxm::KernelType::Cpu4x4_PlainC)
369 reduceGroupEnergySimdBuffers<2>(nbatParams.nenergrp, nbatParams.neg_2log, out);
372 reduceGroupEnergySimdBuffers<4>(nbatParams.nenergrp, nbatParams.neg_2log, out);
375 reduceGroupEnergySimdBuffers<8>(nbatParams.nenergrp, nbatParams.neg_2log, out);
377 default: GMX_RELEASE_ASSERT(false, "Unsupported j-unroll size");
382 wallcycle_sub_stop(wcycle, WallCycleSubCounter::NonbondedKernel);
384 if (stepWork.computeEnergy)
386 reduce_energies_over_lists(nbat, pairlists.ssize(), vVdw, vCoulomb);
390 static void accountFlops(t_nrnb* nrnb,
391 const PairlistSet& pairlistSet,
392 const nonbonded_verlet_t& nbv,
393 const interaction_const_t& ic,
394 const gmx::StepWorkload& stepWork)
396 const bool usingGpuKernels = nbv.useGpu();
398 int enr_nbnxn_kernel_ljc = eNRNB;
399 if (EEL_RF(ic.eeltype) || ic.eeltype == CoulombInteractionType::Cut)
401 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
403 else if ((!usingGpuKernels && nbv.kernelSetup().ewaldExclusionType == Nbnxm::EwaldExclusionType::Analytical)
404 || (usingGpuKernels && Nbnxm::gpu_is_kernel_ewald_analytical(nbv.gpu_nbv)))
406 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
410 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
412 int enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
413 if (stepWork.computeEnergy)
415 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
416 enr_nbnxn_kernel_ljc += 1;
417 enr_nbnxn_kernel_lj += 1;
420 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc, pairlistSet.natpair_ljq_);
421 inc_nrnb(nrnb, enr_nbnxn_kernel_lj, pairlistSet.natpair_lj_);
422 /* The Coulomb-only kernels are offset -eNR_NBNXN_LJ_RF+eNR_NBNXN_RF */
423 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc - eNR_NBNXN_LJ_RF + eNR_NBNXN_RF, pairlistSet.natpair_q_);
425 if (ic.vdw_modifier == InteractionModifiers::ForceSwitch)
427 /* We add up the switch cost separately */
429 eNR_NBNXN_ADD_LJ_FSW + (stepWork.computeEnergy ? 1 : 0),
430 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
432 if (ic.vdw_modifier == InteractionModifiers::PotSwitch)
434 /* We add up the switch cost separately */
436 eNR_NBNXN_ADD_LJ_PSW + (stepWork.computeEnergy ? 1 : 0),
437 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
439 if (ic.vdwtype == VanDerWaalsType::Pme)
441 /* We add up the LJ Ewald cost separately */
443 eNR_NBNXN_ADD_LJ_EWALD + (stepWork.computeEnergy ? 1 : 0),
444 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
448 void nonbonded_verlet_t::dispatchNonbondedKernel(gmx::InteractionLocality iLocality,
449 const interaction_const_t& ic,
450 const gmx::StepWorkload& stepWork,
452 gmx::ArrayRef<const gmx::RVec> shiftvec,
453 gmx::ArrayRef<real> repulsionDispersionSR,
454 gmx::ArrayRef<real> CoulombSR,
457 const PairlistSet& pairlistSet = pairlistSets().pairlistSet(iLocality);
459 switch (kernelSetup().kernelType)
461 case Nbnxm::KernelType::Cpu4x4_PlainC:
462 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
463 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
464 nbnxn_kernel_cpu(pairlistSet,
472 repulsionDispersionSR.data(),
476 case Nbnxm::KernelType::Gpu8x8x8:
477 Nbnxm::gpu_launch_kernel(gpu_nbv, stepWork, iLocality);
480 case Nbnxm::KernelType::Cpu8x8x8_PlainC:
481 nbnxn_kernel_gpu_ref(pairlistSet.gpuList(),
488 nbat->out[0].fshift.data(),
490 repulsionDispersionSR.data());
493 default: GMX_RELEASE_ASSERT(false, "Invalid nonbonded kernel type passed!");
496 accountFlops(nrnb, pairlistSet, *this, ic, stepWork);