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/gmxlib/nonbonded/nb_free_energy.h"
41 #include "gromacs/gmxlib/nonbonded/nb_kernel.h"
42 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
43 #include "gromacs/math/vectypes.h"
44 #include "gromacs/mdlib/enerdata_utils.h"
45 #include "gromacs/mdlib/force.h"
46 #include "gromacs/mdlib/gmx_omp_nthreads.h"
47 #include "gromacs/mdtypes/enerdata.h"
48 #include "gromacs/mdtypes/forceoutput.h"
49 #include "gromacs/mdtypes/forcerec.h"
50 #include "gromacs/mdtypes/inputrec.h"
51 #include "gromacs/mdtypes/interaction_const.h"
52 #include "gromacs/mdtypes/md_enums.h"
53 #include "gromacs/mdtypes/mdatom.h"
54 #include "gromacs/mdtypes/simulation_workload.h"
55 #include "gromacs/nbnxm/gpu_data_mgmt.h"
56 #include "gromacs/nbnxm/nbnxm.h"
57 #include "gromacs/simd/simd.h"
58 #include "gromacs/timing/wallcycle.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_gpu_data_mgmt.h"
65 #include "nbnxm_simd.h"
66 #include "pairlistset.h"
67 #include "pairlistsets.h"
68 #include "kernels_reference/kernel_gpu_ref.h"
69 #define INCLUDE_KERNELFUNCTION_TABLES
70 #include "kernels_reference/kernel_ref.h"
71 #ifdef GMX_NBNXN_SIMD_2XNN
72 # include "kernels_simd_2xmm/kernels.h"
74 #ifdef GMX_NBNXN_SIMD_4XN
75 # include "kernels_simd_4xm/kernels.h"
77 #undef INCLUDE_FUNCTION_TABLES
79 /*! \brief Clears the energy group output buffers
81 * \param[in,out] out nbnxn kernel output struct
83 static void clearGroupEnergies(nbnxn_atomdata_output_t* out)
85 std::fill(out->Vvdw.begin(), out->Vvdw.end(), 0.0_real);
86 std::fill(out->Vc.begin(), out->Vc.end(), 0.0_real);
87 std::fill(out->VSvdw.begin(), out->VSvdw.end(), 0.0_real);
88 std::fill(out->VSc.begin(), out->VSc.end(), 0.0_real);
91 /*! \brief Reduce the group-pair energy buffers produced by a SIMD kernel
92 * to single terms in the output buffers.
94 * The SIMD kernels produce a large number of energy buffer in SIMD registers
95 * to avoid scattered reads and writes.
97 * \tparam unrollj The unroll size for j-particles in the SIMD kernel
98 * \param[in] numGroups The number of energy groups
99 * \param[in] numGroups_2log Log2 of numGroups, rounded up
100 * \param[in,out] out Struct with energy buffers
102 template<int unrollj>
103 static void reduceGroupEnergySimdBuffers(int numGroups, int numGroups_2log, nbnxn_atomdata_output_t* out)
105 const int unrollj_half = unrollj / 2;
106 /* Energies are stored in SIMD registers with size 2^numGroups_2log */
107 const int numGroupsStorage = (1 << numGroups_2log);
109 const real* gmx_restrict vVdwSimd = out->VSvdw.data();
110 const real* gmx_restrict vCoulombSimd = out->VSc.data();
111 real* gmx_restrict vVdw = out->Vvdw.data();
112 real* gmx_restrict vCoulomb = out->Vc.data();
114 /* The size of the SIMD energy group buffer array is:
115 * numGroups*numGroups*numGroupsStorage*unrollj_half*simd_width
117 for (int i = 0; i < numGroups; i++)
119 for (int j1 = 0; j1 < numGroups; j1++)
121 for (int j0 = 0; j0 < numGroups; j0++)
123 int c = ((i * numGroups + j1) * numGroupsStorage + j0) * unrollj_half * unrollj;
124 for (int s = 0; s < unrollj_half; s++)
126 vVdw[i * numGroups + j0] += vVdwSimd[c + 0];
127 vVdw[i * numGroups + j1] += vVdwSimd[c + 1];
128 vCoulomb[i * numGroups + j0] += vCoulombSimd[c + 0];
129 vCoulomb[i * numGroups + j1] += vCoulombSimd[c + 1];
137 /*! \brief Dispatches the non-bonded N versus M atom cluster CPU kernels.
139 * OpenMP parallelization is performed within this function.
140 * Energy reduction, but not force and shift force reduction, is performed
141 * within this function.
143 * \param[in] pairlistSet Pairlists with local or non-local interactions to compute
144 * \param[in] kernelSetup The non-bonded kernel setup
145 * \param[in,out] nbat The atomdata for the interactions
146 * \param[in] ic Non-bonded interaction constants
147 * \param[in] shiftVectors The PBC shift vectors
148 * \param[in] stepWork Flags that tell what to compute
149 * \param[in] clearF Enum that tells if to clear the force output buffer
150 * \param[out] vCoulomb Output buffer for Coulomb energies
151 * \param[out] vVdw Output buffer for Van der Waals energies
152 * \param[in] wcycle Pointer to cycle counting data structure.
154 static void nbnxn_kernel_cpu(const PairlistSet& pairlistSet,
155 const Nbnxm::KernelSetup& kernelSetup,
156 nbnxn_atomdata_t* nbat,
157 const interaction_const_t& ic,
159 const gmx::StepWorkload& stepWork,
163 gmx_wallcycle* wcycle)
167 if (EEL_RF(ic.eeltype) || ic.eeltype == eelCUT)
173 if (kernelSetup.ewaldExclusionType == Nbnxm::EwaldExclusionType::Table)
175 if (ic.rcoulomb == ic.rvdw)
181 coulkt = coulktTAB_TWIN;
186 if (ic.rcoulomb == ic.rvdw)
188 coulkt = coulktEWALD;
192 coulkt = coulktEWALD_TWIN;
197 const nbnxn_atomdata_t::Params& nbatParams = nbat->params();
200 if (ic.vdwtype == evdwCUT)
202 switch (ic.vdw_modifier)
205 case eintmodPOTSHIFT:
206 switch (nbatParams.ljCombinationRule)
208 case LJCombinationRule::Geometric: vdwkt = vdwktLJCUT_COMBGEOM; break;
209 case LJCombinationRule::LorentzBerthelot: vdwkt = vdwktLJCUT_COMBLB; break;
210 case LJCombinationRule::None: vdwkt = vdwktLJCUT_COMBNONE; break;
211 default: GMX_RELEASE_ASSERT(false, "Unknown combination rule");
214 case eintmodFORCESWITCH: vdwkt = vdwktLJFORCESWITCH; break;
215 case eintmodPOTSWITCH: vdwkt = vdwktLJPOTSWITCH; break;
216 default: GMX_RELEASE_ASSERT(false, "Unsupported VdW interaction modifier");
219 else if (ic.vdwtype == evdwPME)
221 if (ic.ljpme_comb_rule == eljpmeGEOM)
223 vdwkt = vdwktLJEWALDCOMBGEOM;
227 vdwkt = vdwktLJEWALDCOMBLB;
228 /* At setup we (should have) selected the C reference kernel */
229 GMX_RELEASE_ASSERT(kernelSetup.kernelType == Nbnxm::KernelType::Cpu4x4_PlainC,
230 "Only the C reference nbnxn SIMD kernel supports LJ-PME with LB "
231 "combination rules");
236 GMX_RELEASE_ASSERT(false, "Unsupported VdW interaction type");
239 gmx::ArrayRef<const NbnxnPairlistCpu> pairlists = pairlistSet.cpuLists();
241 int gmx_unused nthreads = gmx_omp_nthreads_get(emntNonbonded);
242 wallcycle_sub_start(wcycle, ewcsNONBONDED_CLEAR);
243 #pragma omp parallel for schedule(static) num_threads(nthreads)
244 for (gmx::index nb = 0; nb < pairlists.ssize(); nb++)
246 // Presently, the kernels do not call C++ code that can throw,
247 // so no need for a try/catch pair in this OpenMP region.
248 nbnxn_atomdata_output_t* out = &nbat->out[nb];
250 if (clearF == enbvClearFYes)
252 clearForceBuffer(nbat, nb);
254 clear_fshift(out->fshift.data());
259 wallcycle_sub_stop(wcycle, ewcsNONBONDED_CLEAR);
260 wallcycle_sub_start(wcycle, ewcsNONBONDED_KERNEL);
263 // TODO: Change to reference
264 const NbnxnPairlistCpu* pairlist = &pairlists[nb];
266 if (!stepWork.computeEnergy)
268 /* Don't calculate energies */
269 switch (kernelSetup.kernelType)
271 case Nbnxm::KernelType::Cpu4x4_PlainC:
272 nbnxn_kernel_noener_ref[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
274 #ifdef GMX_NBNXN_SIMD_2XNN
275 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
276 nbnxm_kernel_noener_simd_2xmm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
279 #ifdef GMX_NBNXN_SIMD_4XN
280 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
281 nbnxm_kernel_noener_simd_4xm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
284 default: GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
287 else if (out->Vvdw.size() == 1)
289 /* A single energy group (pair) */
293 switch (kernelSetup.kernelType)
295 case Nbnxm::KernelType::Cpu4x4_PlainC:
296 nbnxn_kernel_ener_ref[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
298 #ifdef GMX_NBNXN_SIMD_2XNN
299 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
300 nbnxm_kernel_ener_simd_2xmm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
303 #ifdef GMX_NBNXN_SIMD_4XN
304 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
305 nbnxm_kernel_ener_simd_4xm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
308 default: GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
313 /* Calculate energy group contributions */
314 clearGroupEnergies(out);
318 switch (kernelSetup.kernelType)
320 case Nbnxm::KernelType::Cpu4x4_PlainC:
321 unrollj = c_nbnxnCpuIClusterSize;
322 nbnxn_kernel_energrp_ref[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
324 #ifdef GMX_NBNXN_SIMD_2XNN
325 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
326 unrollj = GMX_SIMD_REAL_WIDTH / 2;
327 nbnxm_kernel_energrp_simd_2xmm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
330 #ifdef GMX_NBNXN_SIMD_4XN
331 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
332 unrollj = GMX_SIMD_REAL_WIDTH;
333 nbnxm_kernel_energrp_simd_4xm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
336 default: GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
339 if (kernelSetup.kernelType != Nbnxm::KernelType::Cpu4x4_PlainC)
344 reduceGroupEnergySimdBuffers<2>(nbatParams.nenergrp, nbatParams.neg_2log, out);
347 reduceGroupEnergySimdBuffers<4>(nbatParams.nenergrp, nbatParams.neg_2log, out);
350 reduceGroupEnergySimdBuffers<8>(nbatParams.nenergrp, nbatParams.neg_2log, out);
352 default: GMX_RELEASE_ASSERT(false, "Unsupported j-unroll size");
357 wallcycle_sub_stop(wcycle, ewcsNONBONDED_KERNEL);
359 if (stepWork.computeEnergy)
361 reduce_energies_over_lists(nbat, pairlists.ssize(), vVdw, vCoulomb);
365 static void accountFlops(t_nrnb* nrnb,
366 const PairlistSet& pairlistSet,
367 const nonbonded_verlet_t& nbv,
368 const interaction_const_t& ic,
369 const gmx::StepWorkload& stepWork)
371 const bool usingGpuKernels = nbv.useGpu();
373 int enr_nbnxn_kernel_ljc;
374 if (EEL_RF(ic.eeltype) || ic.eeltype == eelCUT)
376 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
378 else if ((!usingGpuKernels && nbv.kernelSetup().ewaldExclusionType == Nbnxm::EwaldExclusionType::Analytical)
379 || (usingGpuKernels && Nbnxm::gpu_is_kernel_ewald_analytical(nbv.gpu_nbv)))
381 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
385 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
387 int enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
388 if (stepWork.computeEnergy)
390 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
391 enr_nbnxn_kernel_ljc += 1;
392 enr_nbnxn_kernel_lj += 1;
395 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc, pairlistSet.natpair_ljq_);
396 inc_nrnb(nrnb, enr_nbnxn_kernel_lj, pairlistSet.natpair_lj_);
397 /* The Coulomb-only kernels are offset -eNR_NBNXN_LJ_RF+eNR_NBNXN_RF */
398 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc - eNR_NBNXN_LJ_RF + eNR_NBNXN_RF, pairlistSet.natpair_q_);
400 if (ic.vdw_modifier == eintmodFORCESWITCH)
402 /* We add up the switch cost separately */
404 eNR_NBNXN_ADD_LJ_FSW + (stepWork.computeEnergy ? 1 : 0),
405 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
407 if (ic.vdw_modifier == eintmodPOTSWITCH)
409 /* We add up the switch cost separately */
411 eNR_NBNXN_ADD_LJ_PSW + (stepWork.computeEnergy ? 1 : 0),
412 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
414 if (ic.vdwtype == evdwPME)
416 /* We add up the LJ Ewald cost separately */
418 eNR_NBNXN_ADD_LJ_EWALD + (stepWork.computeEnergy ? 1 : 0),
419 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
423 void nonbonded_verlet_t::dispatchNonbondedKernel(gmx::InteractionLocality iLocality,
424 const interaction_const_t& ic,
425 const gmx::StepWorkload& stepWork,
427 const t_forcerec& fr,
428 gmx_enerdata_t* enerd,
431 const PairlistSet& pairlistSet = pairlistSets().pairlistSet(iLocality);
433 switch (kernelSetup().kernelType)
435 case Nbnxm::KernelType::Cpu4x4_PlainC:
436 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
437 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
438 nbnxn_kernel_cpu(pairlistSet,
445 enerd->grpp.ener[egCOULSR].data(),
446 fr.bBHAM ? enerd->grpp.ener[egBHAMSR].data() : enerd->grpp.ener[egLJSR].data(),
450 case Nbnxm::KernelType::Gpu8x8x8:
451 Nbnxm::gpu_launch_kernel(gpu_nbv, stepWork, iLocality);
454 case Nbnxm::KernelType::Cpu8x8x8_PlainC:
455 nbnxn_kernel_gpu_ref(
456 pairlistSet.gpuList(),
463 nbat->out[0].fshift.data(),
464 enerd->grpp.ener[egCOULSR].data(),
465 fr.bBHAM ? enerd->grpp.ener[egBHAMSR].data() : enerd->grpp.ener[egLJSR].data());
468 default: GMX_RELEASE_ASSERT(false, "Invalid nonbonded kernel type passed!");
471 accountFlops(nrnb, pairlistSet, *this, ic, stepWork);
474 void nonbonded_verlet_t::dispatchFreeEnergyKernel(gmx::InteractionLocality iLocality,
475 const t_forcerec* fr,
477 gmx::ForceWithShiftForces* forceWithShiftForces,
478 const t_mdatoms& mdatoms,
480 gmx::ArrayRef<real const> lambda,
481 gmx_enerdata_t* enerd,
482 const gmx::StepWorkload& stepWork,
485 const auto nbl_fep = pairlistSets().pairlistSet(iLocality).fepLists();
487 /* When the first list is empty, all are empty and there is nothing to do */
488 if (!pairlistSets().params().haveFep || nbl_fep[0]->nrj == 0)
494 /* Add short-range interactions */
495 donb_flags |= GMX_NONBONDED_DO_SR;
497 if (stepWork.computeForces)
499 donb_flags |= GMX_NONBONDED_DO_FORCE;
501 if (stepWork.computeVirial)
503 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
505 if (stepWork.computeEnergy)
507 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
510 nb_kernel_data_t kernel_data;
511 real dvdl_nb[efptNR] = { 0 };
512 kernel_data.flags = donb_flags;
513 kernel_data.lambda = lambda.data();
514 kernel_data.dvdl = dvdl_nb;
516 kernel_data.energygrp_elec = enerd->grpp.ener[egCOULSR].data();
517 kernel_data.energygrp_vdw = enerd->grpp.ener[egLJSR].data();
519 GMX_ASSERT(gmx_omp_nthreads_get(emntNonbonded) == nbl_fep.ssize(),
520 "Number of lists should be same as number of NB threads");
522 wallcycle_sub_start(wcycle_, ewcsNONBONDED_FEP);
523 #pragma omp parallel for schedule(static) num_threads(nbl_fep.ssize())
524 for (gmx::index th = 0; th < nbl_fep.ssize(); th++)
528 gmx_nb_free_energy_kernel(
529 nbl_fep[th].get(), x, forceWithShiftForces, fr, &mdatoms, &kernel_data, nrnb);
531 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
534 if (fepvals->sc_alpha != 0)
536 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
537 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
541 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
542 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
545 /* If we do foreign lambda and we have soft-core interactions
546 * we have to recalculate the (non-linear) energies contributions.
548 if (fepvals->n_lambda > 0 && stepWork.computeDhdl && fepvals->sc_alpha != 0)
551 kernel_data.flags = (donb_flags & ~(GMX_NONBONDED_DO_FORCE | GMX_NONBONDED_DO_SHIFTFORCE))
552 | GMX_NONBONDED_DO_FOREIGNLAMBDA;
553 kernel_data.lambda = lam_i;
554 kernel_data.dvdl = dvdl_nb;
555 kernel_data.energygrp_elec = enerd->foreign_grpp.ener[egCOULSR].data();
556 kernel_data.energygrp_vdw = enerd->foreign_grpp.ener[egLJSR].data();
558 for (gmx::index i = 0; i < 1 + enerd->foreignLambdaTerms.numLambdas(); i++)
560 std::fill(std::begin(dvdl_nb), std::end(dvdl_nb), 0);
561 for (int j = 0; j < efptNR; j++)
563 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i - 1]);
565 reset_foreign_enerdata(enerd);
566 #pragma omp parallel for schedule(static) num_threads(nbl_fep.ssize())
567 for (gmx::index th = 0; th < nbl_fep.ssize(); th++)
571 gmx_nb_free_energy_kernel(
572 nbl_fep[th].get(), x, forceWithShiftForces, fr, &mdatoms, &kernel_data, nrnb);
574 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
577 sum_epot(enerd->foreign_grpp, enerd->foreign_term);
578 enerd->foreignLambdaTerms.accumulate(
579 i, enerd->foreign_term[F_EPOT], dvdl_nb[efptVDW] + dvdl_nb[efptCOUL]);
582 wallcycle_sub_stop(wcycle_, ewcsNONBONDED_FEP);