2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,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.
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.
38 #include "gromacs/gmxlib/nrnb.h"
39 #include "gromacs/gmxlib/nonbonded/nb_free_energy.h"
40 #include "gromacs/gmxlib/nonbonded/nb_kernel.h"
41 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
42 #include "gromacs/math/vectypes.h"
43 #include "gromacs/mdlib/enerdata_utils.h"
44 #include "gromacs/mdlib/force.h"
45 #include "gromacs/mdlib/force_flags.h"
46 #include "gromacs/mdlib/gmx_omp_nthreads.h"
47 #include "gromacs/mdtypes/enerdata.h"
48 #include "gromacs/mdtypes/inputrec.h"
49 #include "gromacs/mdtypes/interaction_const.h"
50 #include "gromacs/mdtypes/md_enums.h"
51 #include "gromacs/mdtypes/mdatom.h"
52 #include "gromacs/nbnxm/gpu_data_mgmt.h"
53 #include "gromacs/nbnxm/nbnxm.h"
54 #include "gromacs/nbnxm/nbnxm_simd.h"
55 #include "gromacs/nbnxm/kernels_reference/kernel_gpu_ref.h"
56 #include "gromacs/simd/simd.h"
57 #include "gromacs/utility/gmxassert.h"
58 #include "gromacs/utility/real.h"
60 #include "kernel_common.h"
61 #include "pairlistset.h"
62 #include "pairlistsets.h"
63 #define INCLUDE_KERNELFUNCTION_TABLES
64 #include "gromacs/nbnxm/kernels_reference/kernel_ref.h"
65 #ifdef GMX_NBNXN_SIMD_2XNN
66 #include "gromacs/nbnxm/kernels_simd_2xmm/kernels.h"
68 #ifdef GMX_NBNXN_SIMD_4XN
69 #include "gromacs/nbnxm/kernels_simd_4xm/kernels.h"
71 #undef INCLUDE_FUNCTION_TABLES
73 /*! \brief Clears the energy group output buffers
75 * \param[in,out] out nbnxn kernel output struct
77 static void clearGroupEnergies(nbnxn_atomdata_output_t *out)
79 std::fill(out->Vvdw.begin(), out->Vvdw.end(), 0.0_real);
80 std::fill(out->Vc.begin(), out->Vc.end(), 0.0_real);
81 std::fill(out->VSvdw.begin(), out->VSvdw.end(), 0.0_real);
82 std::fill(out->VSc.begin(), out->VSc.end(), 0.0_real);
85 /*! \brief Reduce the group-pair energy buffers produced by a SIMD kernel
86 * to single terms in the output buffers.
88 * The SIMD kernels produce a large number of energy buffer in SIMD registers
89 * to avoid scattered reads and writes.
91 * \tparam unrollj The unroll size for j-particles in the SIMD kernel
92 * \param[in] numGroups The number of energy groups
93 * \param[in] numGroups_2log Log2 of numGroups, rounded up
94 * \param[in,out] out Struct with energy buffers
96 template <int unrollj> static void
97 reduceGroupEnergySimdBuffers(int numGroups,
99 nbnxn_atomdata_output_t *out)
101 const int unrollj_half = unrollj/2;
102 /* Energies are stored in SIMD registers with size 2^numGroups_2log */
103 const int numGroupsStorage = (1 << numGroups_2log);
105 const real * gmx_restrict vVdwSimd = out->VSvdw.data();
106 const real * gmx_restrict vCoulombSimd = out->VSc.data();
107 real * gmx_restrict vVdw = out->Vvdw.data();
108 real * gmx_restrict vCoulomb = out->Vc.data();
110 /* The size of the SIMD energy group buffer array is:
111 * numGroups*numGroups*numGroupsStorage*unrollj_half*simd_width
113 for (int i = 0; i < numGroups; i++)
115 for (int j1 = 0; j1 < numGroups; j1++)
117 for (int j0 = 0; j0 < numGroups; j0++)
119 int c = ((i*numGroups + j1)*numGroupsStorage + j0)*unrollj_half*unrollj;
120 for (int s = 0; s < unrollj_half; s++)
122 vVdw [i*numGroups + j0] += vVdwSimd [c + 0];
123 vVdw [i*numGroups + j1] += vVdwSimd [c + 1];
124 vCoulomb[i*numGroups + j0] += vCoulombSimd[c + 0];
125 vCoulomb[i*numGroups + j1] += vCoulombSimd[c + 1];
133 /*! \brief Dispatches the non-bonded N versus M atom cluster CPU kernels.
135 * OpenMP parallelization is performed within this function.
136 * Energy reduction, but not force and shift force reduction, is performed
137 * within this function.
139 * \param[in] pairlistSet Pairlists with local or non-local interactions to compute
140 * \param[in] kernelSetup The non-bonded kernel setup
141 * \param[in,out] nbat The atomdata for the interactions
142 * \param[in] ic Non-bonded interaction constants
143 * \param[in] shiftVectors The PBC shift vectors
144 * \param[in] forceFlags Flags that tell what to compute
145 * \param[in] clearF Enum that tells if to clear the force output buffer
146 * \param[out] vCoulomb Output buffer for Coulomb energies
147 * \param[out] vVdw Output buffer for Van der Waals energies
150 nbnxn_kernel_cpu(const PairlistSet &pairlistSet,
151 const Nbnxm::KernelSetup &kernelSetup,
152 nbnxn_atomdata_t *nbat,
153 const interaction_const_t &ic,
162 if (EEL_RF(ic.eeltype) || ic.eeltype == eelCUT)
168 if (kernelSetup.ewaldExclusionType == Nbnxm::EwaldExclusionType::Table)
170 if (ic.rcoulomb == ic.rvdw)
176 coulkt = coulktTAB_TWIN;
181 if (ic.rcoulomb == ic.rvdw)
183 coulkt = coulktEWALD;
187 coulkt = coulktEWALD_TWIN;
192 const nbnxn_atomdata_t::Params &nbatParams = nbat->params();
195 if (ic.vdwtype == evdwCUT)
197 switch (ic.vdw_modifier)
200 case eintmodPOTSHIFT:
201 switch (nbatParams.comb_rule)
203 case ljcrGEOM: vdwkt = vdwktLJCUT_COMBGEOM; break;
204 case ljcrLB: vdwkt = vdwktLJCUT_COMBLB; break;
205 case ljcrNONE: vdwkt = vdwktLJCUT_COMBNONE; break;
207 GMX_RELEASE_ASSERT(false, "Unknown combination rule");
210 case eintmodFORCESWITCH:
211 vdwkt = vdwktLJFORCESWITCH;
213 case eintmodPOTSWITCH:
214 vdwkt = vdwktLJPOTSWITCH;
217 GMX_RELEASE_ASSERT(false, "Unsupported VdW interaction modifier");
220 else if (ic.vdwtype == evdwPME)
222 if (ic.ljpme_comb_rule == eljpmeGEOM)
224 vdwkt = vdwktLJEWALDCOMBGEOM;
228 vdwkt = vdwktLJEWALDCOMBLB;
229 /* At setup we (should have) selected the C reference kernel */
230 GMX_RELEASE_ASSERT(kernelSetup.kernelType == Nbnxm::KernelType::Cpu4x4_PlainC, "Only the C reference nbnxn SIMD kernel supports LJ-PME with LB combination rules");
235 GMX_RELEASE_ASSERT(false, "Unsupported VdW interaction type");
238 gmx::ArrayRef<const NbnxnPairlistCpu> pairlists = pairlistSet.cpuLists();
240 int gmx_unused nthreads = gmx_omp_nthreads_get(emntNonbonded);
241 #pragma omp parallel for schedule(static) num_threads(nthreads)
242 for (int nb = 0; nb < pairlists.ssize(); nb++)
244 // Presently, the kernels do not call C++ code that can throw,
245 // so no need for a try/catch pair in this OpenMP region.
246 nbnxn_atomdata_output_t *out = &nbat->out[nb];
248 if (clearF == enbvClearFYes)
250 clearForceBuffer(nbat, nb);
252 clear_fshift(out->fshift.data());
255 // TODO: Change to reference
256 const NbnxnPairlistCpu *pairlist = &pairlists[nb];
258 if (!(forceFlags & GMX_FORCE_ENERGY))
260 /* Don't calculate energies */
261 switch (kernelSetup.kernelType)
263 case Nbnxm::KernelType::Cpu4x4_PlainC:
264 nbnxn_kernel_noener_ref[coulkt][vdwkt](pairlist, nbat,
269 #ifdef GMX_NBNXN_SIMD_2XNN
270 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
271 nbnxm_kernel_noener_simd_2xmm[coulkt][vdwkt](pairlist, nbat,
277 #ifdef GMX_NBNXN_SIMD_4XN
278 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
279 nbnxm_kernel_noener_simd_4xm[coulkt][vdwkt](pairlist, nbat,
286 GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
289 else if (out->Vvdw.size() == 1)
291 /* A single energy group (pair) */
295 switch (kernelSetup.kernelType)
297 case Nbnxm::KernelType::Cpu4x4_PlainC:
298 nbnxn_kernel_ener_ref[coulkt][vdwkt](pairlist, nbat,
303 #ifdef GMX_NBNXN_SIMD_2XNN
304 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
305 nbnxm_kernel_ener_simd_2xmm[coulkt][vdwkt](pairlist, nbat,
311 #ifdef GMX_NBNXN_SIMD_4XN
312 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
313 nbnxm_kernel_ener_simd_4xm[coulkt][vdwkt](pairlist, nbat,
320 GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
325 /* Calculate energy group contributions */
326 clearGroupEnergies(out);
330 switch (kernelSetup.kernelType)
332 case Nbnxm::KernelType::Cpu4x4_PlainC:
333 unrollj = c_nbnxnCpuIClusterSize;
334 nbnxn_kernel_energrp_ref[coulkt][vdwkt](pairlist, nbat,
339 #ifdef GMX_NBNXN_SIMD_2XNN
340 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
341 unrollj = GMX_SIMD_REAL_WIDTH/2;
342 nbnxm_kernel_energrp_simd_2xmm[coulkt][vdwkt](pairlist, nbat,
348 #ifdef GMX_NBNXN_SIMD_4XN
349 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
350 unrollj = GMX_SIMD_REAL_WIDTH;
351 nbnxm_kernel_energrp_simd_4xm[coulkt][vdwkt](pairlist, nbat,
358 GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
361 if (kernelSetup.kernelType != Nbnxm::KernelType::Cpu4x4_PlainC)
366 reduceGroupEnergySimdBuffers<2>(nbatParams.nenergrp,
371 reduceGroupEnergySimdBuffers<4>(nbatParams.nenergrp,
376 reduceGroupEnergySimdBuffers<8>(nbatParams.nenergrp,
381 GMX_RELEASE_ASSERT(false, "Unsupported j-unroll size");
387 if (forceFlags & GMX_FORCE_ENERGY)
389 reduce_energies_over_lists(nbat, pairlists.ssize(), vVdw, vCoulomb);
393 static void accountFlops(t_nrnb *nrnb,
394 const PairlistSet &pairlistSet,
395 const nonbonded_verlet_t &nbv,
396 const interaction_const_t &ic,
397 const int forceFlags)
399 const bool usingGpuKernels = nbv.useGpu();
401 int enr_nbnxn_kernel_ljc;
402 if (EEL_RF(ic.eeltype) || ic.eeltype == eelCUT)
404 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
406 else if ((!usingGpuKernels && nbv.kernelSetup().ewaldExclusionType == Nbnxm::EwaldExclusionType::Analytical) ||
407 (usingGpuKernels && Nbnxm::gpu_is_kernel_ewald_analytical(nbv.gpu_nbv)))
409 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
413 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
415 int enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
416 if (forceFlags & GMX_FORCE_ENERGY)
418 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
419 enr_nbnxn_kernel_ljc += 1;
420 enr_nbnxn_kernel_lj += 1;
423 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc,
424 pairlistSet.natpair_ljq_);
425 inc_nrnb(nrnb, enr_nbnxn_kernel_lj,
426 pairlistSet.natpair_lj_);
427 /* The Coulomb-only kernels are offset -eNR_NBNXN_LJ_RF+eNR_NBNXN_RF */
428 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
429 pairlistSet.natpair_q_);
431 const bool calcEnergy = ((forceFlags & GMX_FORCE_ENERGY) != 0);
432 if (ic.vdw_modifier == eintmodFORCESWITCH)
434 /* We add up the switch cost separately */
435 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_FSW + (calcEnergy ? 1 : 0),
436 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
438 if (ic.vdw_modifier == eintmodPOTSWITCH)
440 /* We add up the switch cost separately */
441 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_PSW + (calcEnergy ? 1 : 0),
442 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
444 if (ic.vdwtype == evdwPME)
446 /* We add up the LJ Ewald cost separately */
447 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_EWALD + (calcEnergy ? 1 : 0),
448 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
453 nonbonded_verlet_t::dispatchNonbondedKernel(Nbnxm::InteractionLocality iLocality,
454 const interaction_const_t &ic,
457 const t_forcerec &fr,
458 gmx_enerdata_t *enerd,
461 const PairlistSet &pairlistSet = pairlistSets().pairlistSet(iLocality);
463 switch (kernelSetup().kernelType)
465 case Nbnxm::KernelType::Cpu4x4_PlainC:
466 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
467 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
468 nbnxn_kernel_cpu(pairlistSet,
475 enerd->grpp.ener[egCOULSR].data(),
477 enerd->grpp.ener[egBHAMSR].data() :
478 enerd->grpp.ener[egLJSR].data());
481 case Nbnxm::KernelType::Gpu8x8x8:
482 Nbnxm::gpu_launch_kernel(gpu_nbv, forceFlags, iLocality);
485 case Nbnxm::KernelType::Cpu8x8x8_PlainC:
486 nbnxn_kernel_gpu_ref(pairlistSet.gpuList(),
492 nbat->out[0].fshift.data(),
493 enerd->grpp.ener[egCOULSR].data(),
495 enerd->grpp.ener[egBHAMSR].data() :
496 enerd->grpp.ener[egLJSR].data());
500 GMX_RELEASE_ASSERT(false, "Invalid nonbonded kernel type passed!");
504 accountFlops(nrnb, pairlistSet, *this, ic, forceFlags);
508 nonbonded_verlet_t::dispatchFreeEnergyKernel(Nbnxm::InteractionLocality iLocality,
512 const t_mdatoms &mdatoms,
515 gmx_enerdata_t *enerd,
516 const int forceFlags,
519 const auto nbl_fep = pairlistSets().pairlistSet(iLocality).fepLists();
521 /* When the first list is empty, all are empty and there is nothing to do */
522 if (!pairlistSets().params().haveFep || nbl_fep[0]->nrj == 0)
528 /* Add short-range interactions */
529 donb_flags |= GMX_NONBONDED_DO_SR;
531 /* Currently all group scheme kernels always calculate (shift-)forces */
532 if (forceFlags & GMX_FORCE_FORCES)
534 donb_flags |= GMX_NONBONDED_DO_FORCE;
536 if (forceFlags & GMX_FORCE_VIRIAL)
538 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
540 if (forceFlags & GMX_FORCE_ENERGY)
542 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
545 nb_kernel_data_t kernel_data;
546 real dvdl_nb[efptNR] = { 0 };
547 kernel_data.flags = donb_flags;
548 kernel_data.lambda = lambda;
549 kernel_data.dvdl = dvdl_nb;
551 kernel_data.energygrp_elec = enerd->grpp.ener[egCOULSR].data();
552 kernel_data.energygrp_vdw = enerd->grpp.ener[egLJSR].data();
554 GMX_ASSERT(gmx_omp_nthreads_get(emntNonbonded) == nbl_fep.ssize(), "Number of lists should be same as number of NB threads");
556 #pragma omp parallel for schedule(static) num_threads(nbl_fep.ssize())
557 for (int th = 0; th < nbl_fep.ssize(); th++)
561 gmx_nb_free_energy_kernel(nbl_fep[th].get(),
562 x, f, fr, &mdatoms, &kernel_data, nrnb);
564 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
567 if (fepvals->sc_alpha != 0)
569 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
570 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
574 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
575 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
578 /* If we do foreign lambda and we have soft-core interactions
579 * we have to recalculate the (non-linear) energies contributions.
581 if (fepvals->n_lambda > 0 && (forceFlags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
584 kernel_data.flags = (donb_flags & ~(GMX_NONBONDED_DO_FORCE | GMX_NONBONDED_DO_SHIFTFORCE)) | GMX_NONBONDED_DO_FOREIGNLAMBDA;
585 kernel_data.lambda = lam_i;
586 kernel_data.energygrp_elec = enerd->foreign_grpp.ener[egCOULSR].data();
587 kernel_data.energygrp_vdw = enerd->foreign_grpp.ener[egLJSR].data();
588 /* Note that we add to kernel_data.dvdl, but ignore the result */
590 for (size_t i = 0; i < enerd->enerpart_lambda.size(); i++)
592 for (int j = 0; j < efptNR; j++)
594 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
596 reset_foreign_enerdata(enerd);
597 #pragma omp parallel for schedule(static) num_threads(nbl_fep.ssize())
598 for (int th = 0; th < nbl_fep.ssize(); th++)
602 gmx_nb_free_energy_kernel(nbl_fep[th].get(),
603 x, f, fr, &mdatoms, &kernel_data, nrnb);
605 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
608 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
609 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];