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/force.h"
44 #include "gromacs/mdlib/force_flags.h"
45 #include "gromacs/mdlib/gmx_omp_nthreads.h"
46 #include "gromacs/mdtypes/enerdata.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/nbnxm/gpu_data_mgmt.h"
52 #include "gromacs/nbnxm/nbnxm.h"
53 #include "gromacs/nbnxm/nbnxm_simd.h"
54 #include "gromacs/nbnxm/kernels_reference/kernel_gpu_ref.h"
55 #include "gromacs/simd/simd.h"
56 #include "gromacs/utility/gmxassert.h"
57 #include "gromacs/utility/real.h"
59 #include "kernel_common.h"
60 #define INCLUDE_KERNELFUNCTION_TABLES
61 #include "gromacs/nbnxm/kernels_reference/kernel_ref.h"
62 #ifdef GMX_NBNXN_SIMD_2XNN
63 #include "gromacs/nbnxm/kernels_simd_2xmm/kernels.h"
65 #ifdef GMX_NBNXN_SIMD_4XN
66 #include "gromacs/nbnxm/kernels_simd_4xm/kernels.h"
68 #undef INCLUDE_FUNCTION_TABLES
70 /*! \brief Clears the energy group output buffers
72 * \param[in,out] out nbnxn kernel output struct
74 static void clearGroupEnergies(nbnxn_atomdata_output_t *out)
76 std::fill(out->Vvdw.begin(), out->Vvdw.end(), 0.0_real);
77 std::fill(out->Vc.begin(), out->Vc.end(), 0.0_real);
78 std::fill(out->VSvdw.begin(), out->VSvdw.end(), 0.0_real);
79 std::fill(out->VSc.begin(), out->VSc.end(), 0.0_real);
82 /*! \brief Reduce the group-pair energy buffers produced by a SIMD kernel
83 * to single terms in the output buffers.
85 * The SIMD kernels produce a large number of energy buffer in SIMD registers
86 * to avoid scattered reads and writes.
88 * \tparam unrollj The unroll size for j-particles in the SIMD kernel
89 * \param[in] numGroups The number of energy groups
90 * \param[in] numGroups_2log Log2 of numGroups, rounded up
91 * \param[in,out] out Struct with energy buffers
93 template <int unrollj> static void
94 reduceGroupEnergySimdBuffers(int numGroups,
96 nbnxn_atomdata_output_t *out)
98 const int unrollj_half = unrollj/2;
99 /* Energies are stored in SIMD registers with size 2^numGroups_2log */
100 const int numGroupsStorage = (1 << numGroups_2log);
102 const real * gmx_restrict vVdwSimd = out->VSvdw.data();
103 const real * gmx_restrict vCoulombSimd = out->VSc.data();
104 real * gmx_restrict vVdw = out->Vvdw.data();
105 real * gmx_restrict vCoulomb = out->Vc.data();
107 /* The size of the SIMD energy group buffer array is:
108 * numGroups*numGroups*numGroupsStorage*unrollj_half*simd_width
110 for (int i = 0; i < numGroups; i++)
112 for (int j1 = 0; j1 < numGroups; j1++)
114 for (int j0 = 0; j0 < numGroups; j0++)
116 int c = ((i*numGroups + j1)*numGroupsStorage + j0)*unrollj_half*unrollj;
117 for (int s = 0; s < unrollj_half; s++)
119 vVdw [i*numGroups + j0] += vVdwSimd [c + 0];
120 vVdw [i*numGroups + j1] += vVdwSimd [c + 1];
121 vCoulomb[i*numGroups + j0] += vCoulombSimd[c + 0];
122 vCoulomb[i*numGroups + j1] += vCoulombSimd[c + 1];
130 /*! \brief Dispatches the non-bonded N versus M atom cluster CPU kernels.
132 * OpenMP parallelization is performed within this function.
133 * Energy reduction, but not force and shift force reduction, is performed
134 * within this function.
136 * \param[in] pairlistSet Pairlists with local or non-local interactions to compute
137 * \param[in] kernelSetup The non-bonded kernel setup
138 * \param[in,out] nbat The atomdata for the interactions
139 * \param[in] ic Non-bonded interaction constants
140 * \param[in] shiftVectors The PBC shift vectors
141 * \param[in] forceFlags Flags that tell what to compute
142 * \param[in] clearF Enum that tells if to clear the force output buffer
143 * \param[out] fshift Shift force output buffer
144 * \param[out] vCoulomb Output buffer for Coulomb energies
145 * \param[out] vVdw Output buffer for Van der Waals energies
148 nbnxn_kernel_cpu(const nbnxn_pairlist_set_t &pairlistSet,
149 const Nbnxm::KernelSetup &kernelSetup,
150 nbnxn_atomdata_t *nbat,
151 const interaction_const_t &ic,
161 if (EEL_RF(ic.eeltype) || ic.eeltype == eelCUT)
167 if (kernelSetup.ewaldExclusionType == Nbnxm::EwaldExclusionType::Table)
169 if (ic.rcoulomb == ic.rvdw)
175 coulkt = coulktTAB_TWIN;
180 if (ic.rcoulomb == ic.rvdw)
182 coulkt = coulktEWALD;
186 coulkt = coulktEWALD_TWIN;
191 const nbnxn_atomdata_t::Params &nbatParams = nbat->params();
194 if (ic.vdwtype == evdwCUT)
196 switch (ic.vdw_modifier)
199 case eintmodPOTSHIFT:
200 switch (nbatParams.comb_rule)
202 case ljcrGEOM: vdwkt = vdwktLJCUT_COMBGEOM; break;
203 case ljcrLB: vdwkt = vdwktLJCUT_COMBLB; break;
204 case ljcrNONE: vdwkt = vdwktLJCUT_COMBNONE; break;
206 GMX_RELEASE_ASSERT(false, "Unknown combination rule");
209 case eintmodFORCESWITCH:
210 vdwkt = vdwktLJFORCESWITCH;
212 case eintmodPOTSWITCH:
213 vdwkt = vdwktLJPOTSWITCH;
216 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, "Only the C reference nbnxn SIMD kernel supports LJ-PME with LB combination rules");
234 GMX_RELEASE_ASSERT(false, "Unsupported VdW interaction type");
237 int nnbl = pairlistSet.nnbl;
238 NbnxnPairlistCpu * const * nbl = pairlistSet.nbl;
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 < nnbl; 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 clear_f(nbat, nb, out->f.data());
254 if ((forceFlags & GMX_FORCE_VIRIAL) && nnbl == 1)
260 fshift_p = out->fshift.data();
262 if (clearF == enbvClearFYes)
264 clear_fshift(fshift_p);
268 if (!(forceFlags & GMX_FORCE_ENERGY))
270 /* Don't calculate energies */
271 switch (kernelSetup.kernelType)
273 case Nbnxm::KernelType::Cpu4x4_PlainC:
274 nbnxn_kernel_noener_ref[coulkt][vdwkt](nbl[nb], nbat,
280 #ifdef GMX_NBNXN_SIMD_2XNN
281 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
282 nbnxm_kernel_noener_simd_2xmm[coulkt][vdwkt](nbl[nb], nbat,
289 #ifdef GMX_NBNXN_SIMD_4XN
290 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
291 nbnxm_kernel_noener_simd_4xm[coulkt][vdwkt](nbl[nb], nbat,
299 GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
302 else if (out->Vvdw.size() == 1)
304 /* A single energy group (pair) */
308 switch (kernelSetup.kernelType)
310 case Nbnxm::KernelType::Cpu4x4_PlainC:
311 nbnxn_kernel_ener_ref[coulkt][vdwkt](nbl[nb], nbat,
319 #ifdef GMX_NBNXN_SIMD_2XNN
320 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
321 nbnxm_kernel_ener_simd_2xmm[coulkt][vdwkt](nbl[nb], nbat,
330 #ifdef GMX_NBNXN_SIMD_4XN
331 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
332 nbnxm_kernel_ener_simd_4xm[coulkt][vdwkt](nbl[nb], nbat,
342 GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
347 /* Calculate energy group contributions */
348 clearGroupEnergies(out);
352 switch (kernelSetup.kernelType)
354 case Nbnxm::KernelType::Cpu4x4_PlainC:
355 unrollj = c_nbnxnCpuIClusterSize;
356 nbnxn_kernel_energrp_ref[coulkt][vdwkt](nbl[nb], nbat,
364 #ifdef GMX_NBNXN_SIMD_2XNN
365 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
366 unrollj = GMX_SIMD_REAL_WIDTH/2;
367 nbnxm_kernel_energrp_simd_2xmm[coulkt][vdwkt](nbl[nb], nbat,
376 #ifdef GMX_NBNXN_SIMD_4XN
377 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
378 unrollj = GMX_SIMD_REAL_WIDTH;
379 nbnxm_kernel_energrp_simd_4xm[coulkt][vdwkt](nbl[nb], nbat,
389 GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
392 if (kernelSetup.kernelType != Nbnxm::KernelType::Cpu4x4_PlainC)
397 reduceGroupEnergySimdBuffers<2>(nbatParams.nenergrp,
402 reduceGroupEnergySimdBuffers<4>(nbatParams.nenergrp,
407 reduceGroupEnergySimdBuffers<8>(nbatParams.nenergrp,
412 GMX_RELEASE_ASSERT(false, "Unsupported j-unroll size");
418 if (forceFlags & GMX_FORCE_ENERGY)
420 reduce_energies_over_lists(nbat, nnbl, vVdw, vCoulomb);
424 static void accountFlops(t_nrnb *nrnb,
425 const nonbonded_verlet_t &nbv,
426 const Nbnxm::InteractionLocality iLocality,
427 const interaction_const_t &ic,
428 const int forceFlags)
430 const nbnxn_pairlist_set_t &pairlistSet = nbv.pairlistSet(iLocality);
431 const bool usingGpuKernels = nbv.useGpu();
433 int enr_nbnxn_kernel_ljc;
434 if (EEL_RF(ic.eeltype) || ic.eeltype == eelCUT)
436 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
438 else if ((!usingGpuKernels && nbv.kernelSetup().ewaldExclusionType == Nbnxm::EwaldExclusionType::Analytical) ||
439 (usingGpuKernels && Nbnxm::gpu_is_kernel_ewald_analytical(nbv.gpu_nbv)))
441 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
445 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
447 int enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
448 if (forceFlags & GMX_FORCE_ENERGY)
450 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
451 enr_nbnxn_kernel_ljc += 1;
452 enr_nbnxn_kernel_lj += 1;
455 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc,
456 pairlistSet.natpair_ljq);
457 inc_nrnb(nrnb, enr_nbnxn_kernel_lj,
458 pairlistSet.natpair_lj);
459 /* The Coulomb-only kernels are offset -eNR_NBNXN_LJ_RF+eNR_NBNXN_RF */
460 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
461 pairlistSet.natpair_q);
463 const bool calcEnergy = ((forceFlags & GMX_FORCE_ENERGY) != 0);
464 if (ic.vdw_modifier == eintmodFORCESWITCH)
466 /* We add up the switch cost separately */
467 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_FSW + (calcEnergy ? 1 : 0),
468 pairlistSet.natpair_ljq + pairlistSet.natpair_lj);
470 if (ic.vdw_modifier == eintmodPOTSWITCH)
472 /* We add up the switch cost separately */
473 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_PSW + (calcEnergy ? 1 : 0),
474 pairlistSet.natpair_ljq + pairlistSet.natpair_lj);
476 if (ic.vdwtype == evdwPME)
478 /* We add up the LJ Ewald cost separately */
479 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_EWALD + (calcEnergy ? 1 : 0),
480 pairlistSet.natpair_ljq + pairlistSet.natpair_lj);
484 void NbnxnDispatchKernel(nonbonded_verlet_t *nbv,
485 Nbnxm::InteractionLocality iLocality,
486 const interaction_const_t &ic,
490 gmx_enerdata_t *enerd,
493 const nbnxn_pairlist_set_t &pairlistSet = nbv->pairlistSet(iLocality);
495 switch (nbv->kernelSetup().kernelType)
497 case Nbnxm::KernelType::Cpu4x4_PlainC:
498 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
499 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
500 nbnxn_kernel_cpu(pairlistSet,
508 enerd->grpp.ener[egCOULSR],
510 enerd->grpp.ener[egBHAMSR] :
511 enerd->grpp.ener[egLJSR]);
514 case Nbnxm::KernelType::Gpu8x8x8:
515 Nbnxm::gpu_launch_kernel(nbv->gpu_nbv, forceFlags, iLocality);
518 case Nbnxm::KernelType::Cpu8x8x8_PlainC:
519 nbnxn_kernel_gpu_ref(pairlistSet.nblGpu[0],
526 enerd->grpp.ener[egCOULSR],
528 enerd->grpp.ener[egBHAMSR] :
529 enerd->grpp.ener[egLJSR]);
533 GMX_RELEASE_ASSERT(false, "Invalid nonbonded kernel type passed!");
537 accountFlops(nrnb, *nbv, iLocality, ic, forceFlags);
541 nonbonded_verlet_t::dispatchFreeEnergyKernel(Nbnxm::InteractionLocality iLocality,
545 const t_mdatoms &mdatoms,
548 gmx_enerdata_t *enerd,
549 const int forceFlags,
552 const gmx::ArrayRef<t_nblist const * const > nbl_fep = pairlistSet(iLocality).nbl_fep;
554 /* When the first list is empty, all are empty and there is nothing to do */
555 if (nbl_fep[0]->nrj == 0)
561 /* Add short-range interactions */
562 donb_flags |= GMX_NONBONDED_DO_SR;
564 /* Currently all group scheme kernels always calculate (shift-)forces */
565 if (forceFlags & GMX_FORCE_FORCES)
567 donb_flags |= GMX_NONBONDED_DO_FORCE;
569 if (forceFlags & GMX_FORCE_VIRIAL)
571 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
573 if (forceFlags & GMX_FORCE_ENERGY)
575 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
578 nb_kernel_data_t kernel_data;
579 real dvdl_nb[efptNR] = { 0 };
580 kernel_data.flags = donb_flags;
581 kernel_data.lambda = lambda;
582 kernel_data.dvdl = dvdl_nb;
584 kernel_data.energygrp_elec = enerd->grpp.ener[egCOULSR];
585 kernel_data.energygrp_vdw = enerd->grpp.ener[egLJSR];
587 GMX_ASSERT(gmx_omp_nthreads_get(emntNonbonded) == nbl_fep.ssize(), "Number of lists should be same as number of NB threads");
589 #pragma omp parallel for schedule(static) num_threads(nbl_fep.ssize())
590 for (int th = 0; th < nbl_fep.ssize(); th++)
594 gmx_nb_free_energy_kernel(nbl_fep[th],
595 x, f, fr, &mdatoms, &kernel_data, nrnb);
597 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
600 if (fepvals->sc_alpha != 0)
602 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
603 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
607 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
608 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
611 /* If we do foreign lambda and we have soft-core interactions
612 * we have to recalculate the (non-linear) energies contributions.
614 if (fepvals->n_lambda > 0 && (forceFlags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
617 kernel_data.flags = (donb_flags & ~(GMX_NONBONDED_DO_FORCE | GMX_NONBONDED_DO_SHIFTFORCE)) | GMX_NONBONDED_DO_FOREIGNLAMBDA;
618 kernel_data.lambda = lam_i;
619 kernel_data.energygrp_elec = enerd->foreign_grpp.ener[egCOULSR];
620 kernel_data.energygrp_vdw = enerd->foreign_grpp.ener[egLJSR];
621 /* Note that we add to kernel_data.dvdl, but ignore the result */
623 for (int i = 0; i < enerd->n_lambda; i++)
625 for (int j = 0; j < efptNR; j++)
627 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
629 reset_foreign_enerdata(enerd);
630 #pragma omp parallel for schedule(static) num_threads(nbl_fep.ssize())
631 for (int th = 0; th < nbl_fep.ssize(); th++)
635 gmx_nb_free_energy_kernel(nbl_fep[th],
636 x, f, fr, &mdatoms, &kernel_data, nrnb);
638 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
641 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
642 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];