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 #include "pairlistset.h"
61 #include "pairlistsets.h"
62 #define INCLUDE_KERNELFUNCTION_TABLES
63 #include "gromacs/nbnxm/kernels_reference/kernel_ref.h"
64 #ifdef GMX_NBNXN_SIMD_2XNN
65 #include "gromacs/nbnxm/kernels_simd_2xmm/kernels.h"
67 #ifdef GMX_NBNXN_SIMD_4XN
68 #include "gromacs/nbnxm/kernels_simd_4xm/kernels.h"
70 #undef INCLUDE_FUNCTION_TABLES
72 /*! \brief Clears the energy group output buffers
74 * \param[in,out] out nbnxn kernel output struct
76 static void clearGroupEnergies(nbnxn_atomdata_output_t *out)
78 std::fill(out->Vvdw.begin(), out->Vvdw.end(), 0.0_real);
79 std::fill(out->Vc.begin(), out->Vc.end(), 0.0_real);
80 std::fill(out->VSvdw.begin(), out->VSvdw.end(), 0.0_real);
81 std::fill(out->VSc.begin(), out->VSc.end(), 0.0_real);
84 /*! \brief Reduce the group-pair energy buffers produced by a SIMD kernel
85 * to single terms in the output buffers.
87 * The SIMD kernels produce a large number of energy buffer in SIMD registers
88 * to avoid scattered reads and writes.
90 * \tparam unrollj The unroll size for j-particles in the SIMD kernel
91 * \param[in] numGroups The number of energy groups
92 * \param[in] numGroups_2log Log2 of numGroups, rounded up
93 * \param[in,out] out Struct with energy buffers
95 template <int unrollj> static void
96 reduceGroupEnergySimdBuffers(int numGroups,
98 nbnxn_atomdata_output_t *out)
100 const int unrollj_half = unrollj/2;
101 /* Energies are stored in SIMD registers with size 2^numGroups_2log */
102 const int numGroupsStorage = (1 << numGroups_2log);
104 const real * gmx_restrict vVdwSimd = out->VSvdw.data();
105 const real * gmx_restrict vCoulombSimd = out->VSc.data();
106 real * gmx_restrict vVdw = out->Vvdw.data();
107 real * gmx_restrict vCoulomb = out->Vc.data();
109 /* The size of the SIMD energy group buffer array is:
110 * numGroups*numGroups*numGroupsStorage*unrollj_half*simd_width
112 for (int i = 0; i < numGroups; i++)
114 for (int j1 = 0; j1 < numGroups; j1++)
116 for (int j0 = 0; j0 < numGroups; j0++)
118 int c = ((i*numGroups + j1)*numGroupsStorage + j0)*unrollj_half*unrollj;
119 for (int s = 0; s < unrollj_half; s++)
121 vVdw [i*numGroups + j0] += vVdwSimd [c + 0];
122 vVdw [i*numGroups + j1] += vVdwSimd [c + 1];
123 vCoulomb[i*numGroups + j0] += vCoulombSimd[c + 0];
124 vCoulomb[i*numGroups + j1] += vCoulombSimd[c + 1];
132 /*! \brief Dispatches the non-bonded N versus M atom cluster CPU kernels.
134 * OpenMP parallelization is performed within this function.
135 * Energy reduction, but not force and shift force reduction, is performed
136 * within this function.
138 * \param[in] pairlistSet Pairlists with local or non-local interactions to compute
139 * \param[in] kernelSetup The non-bonded kernel setup
140 * \param[in,out] nbat The atomdata for the interactions
141 * \param[in] ic Non-bonded interaction constants
142 * \param[in] shiftVectors The PBC shift vectors
143 * \param[in] forceFlags Flags that tell what to compute
144 * \param[in] clearF Enum that tells if to clear the force output buffer
145 * \param[out] fshift Shift 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,
163 if (EEL_RF(ic.eeltype) || ic.eeltype == eelCUT)
169 if (kernelSetup.ewaldExclusionType == Nbnxm::EwaldExclusionType::Table)
171 if (ic.rcoulomb == ic.rvdw)
177 coulkt = coulktTAB_TWIN;
182 if (ic.rcoulomb == ic.rvdw)
184 coulkt = coulktEWALD;
188 coulkt = coulktEWALD_TWIN;
193 const nbnxn_atomdata_t::Params &nbatParams = nbat->params();
196 if (ic.vdwtype == evdwCUT)
198 switch (ic.vdw_modifier)
201 case eintmodPOTSHIFT:
202 switch (nbatParams.comb_rule)
204 case ljcrGEOM: vdwkt = vdwktLJCUT_COMBGEOM; break;
205 case ljcrLB: vdwkt = vdwktLJCUT_COMBLB; break;
206 case ljcrNONE: vdwkt = vdwktLJCUT_COMBNONE; break;
208 GMX_RELEASE_ASSERT(false, "Unknown combination rule");
211 case eintmodFORCESWITCH:
212 vdwkt = vdwktLJFORCESWITCH;
214 case eintmodPOTSWITCH:
215 vdwkt = vdwktLJPOTSWITCH;
218 GMX_RELEASE_ASSERT(false, "Unsupported VdW interaction modifier");
221 else if (ic.vdwtype == evdwPME)
223 if (ic.ljpme_comb_rule == eljpmeGEOM)
225 vdwkt = vdwktLJEWALDCOMBGEOM;
229 vdwkt = vdwktLJEWALDCOMBLB;
230 /* At setup we (should have) selected the C reference kernel */
231 GMX_RELEASE_ASSERT(kernelSetup.kernelType == Nbnxm::KernelType::Cpu4x4_PlainC, "Only the C reference nbnxn SIMD kernel supports LJ-PME with LB 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 #pragma omp parallel for schedule(static) num_threads(nthreads)
243 for (int nb = 0; nb < pairlists.ssize(); nb++)
245 // Presently, the kernels do not call C++ code that can throw,
246 // so no need for a try/catch pair in this OpenMP region.
247 nbnxn_atomdata_output_t *out = &nbat->out[nb];
249 if (clearF == enbvClearFYes)
251 clear_f(nbat, nb, out->f.data());
255 if ((forceFlags & GMX_FORCE_VIRIAL) && pairlists.ssize() == 1)
261 fshift_p = out->fshift.data();
263 if (clearF == enbvClearFYes)
265 clear_fshift(fshift_p);
269 // TODO: Change to reference
270 const NbnxnPairlistCpu *pairlist = &pairlists[nb];
272 if (!(forceFlags & GMX_FORCE_ENERGY))
274 /* Don't calculate energies */
275 switch (kernelSetup.kernelType)
277 case Nbnxm::KernelType::Cpu4x4_PlainC:
278 nbnxn_kernel_noener_ref[coulkt][vdwkt](pairlist, nbat,
284 #ifdef GMX_NBNXN_SIMD_2XNN
285 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
286 nbnxm_kernel_noener_simd_2xmm[coulkt][vdwkt](pairlist, nbat,
293 #ifdef GMX_NBNXN_SIMD_4XN
294 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
295 nbnxm_kernel_noener_simd_4xm[coulkt][vdwkt](pairlist, nbat,
303 GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
306 else if (out->Vvdw.size() == 1)
308 /* A single energy group (pair) */
312 switch (kernelSetup.kernelType)
314 case Nbnxm::KernelType::Cpu4x4_PlainC:
315 nbnxn_kernel_ener_ref[coulkt][vdwkt](pairlist, nbat,
323 #ifdef GMX_NBNXN_SIMD_2XNN
324 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
325 nbnxm_kernel_ener_simd_2xmm[coulkt][vdwkt](pairlist, nbat,
334 #ifdef GMX_NBNXN_SIMD_4XN
335 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
336 nbnxm_kernel_ener_simd_4xm[coulkt][vdwkt](pairlist, nbat,
346 GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
351 /* Calculate energy group contributions */
352 clearGroupEnergies(out);
356 switch (kernelSetup.kernelType)
358 case Nbnxm::KernelType::Cpu4x4_PlainC:
359 unrollj = c_nbnxnCpuIClusterSize;
360 nbnxn_kernel_energrp_ref[coulkt][vdwkt](pairlist, nbat,
368 #ifdef GMX_NBNXN_SIMD_2XNN
369 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
370 unrollj = GMX_SIMD_REAL_WIDTH/2;
371 nbnxm_kernel_energrp_simd_2xmm[coulkt][vdwkt](pairlist, nbat,
380 #ifdef GMX_NBNXN_SIMD_4XN
381 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
382 unrollj = GMX_SIMD_REAL_WIDTH;
383 nbnxm_kernel_energrp_simd_4xm[coulkt][vdwkt](pairlist, nbat,
393 GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
396 if (kernelSetup.kernelType != Nbnxm::KernelType::Cpu4x4_PlainC)
401 reduceGroupEnergySimdBuffers<2>(nbatParams.nenergrp,
406 reduceGroupEnergySimdBuffers<4>(nbatParams.nenergrp,
411 reduceGroupEnergySimdBuffers<8>(nbatParams.nenergrp,
416 GMX_RELEASE_ASSERT(false, "Unsupported j-unroll size");
422 if (forceFlags & GMX_FORCE_ENERGY)
424 reduce_energies_over_lists(nbat, pairlists.ssize(), vVdw, vCoulomb);
428 static void accountFlops(t_nrnb *nrnb,
429 const PairlistSet &pairlistSet,
430 const nonbonded_verlet_t &nbv,
431 const interaction_const_t &ic,
432 const int forceFlags)
434 const bool usingGpuKernels = nbv.useGpu();
436 int enr_nbnxn_kernel_ljc;
437 if (EEL_RF(ic.eeltype) || ic.eeltype == eelCUT)
439 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
441 else if ((!usingGpuKernels && nbv.kernelSetup().ewaldExclusionType == Nbnxm::EwaldExclusionType::Analytical) ||
442 (usingGpuKernels && Nbnxm::gpu_is_kernel_ewald_analytical(nbv.gpu_nbv)))
444 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
448 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
450 int enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
451 if (forceFlags & GMX_FORCE_ENERGY)
453 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
454 enr_nbnxn_kernel_ljc += 1;
455 enr_nbnxn_kernel_lj += 1;
458 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc,
459 pairlistSet.natpair_ljq_);
460 inc_nrnb(nrnb, enr_nbnxn_kernel_lj,
461 pairlistSet.natpair_lj_);
462 /* The Coulomb-only kernels are offset -eNR_NBNXN_LJ_RF+eNR_NBNXN_RF */
463 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
464 pairlistSet.natpair_q_);
466 const bool calcEnergy = ((forceFlags & GMX_FORCE_ENERGY) != 0);
467 if (ic.vdw_modifier == eintmodFORCESWITCH)
469 /* We add up the switch cost separately */
470 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_FSW + (calcEnergy ? 1 : 0),
471 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
473 if (ic.vdw_modifier == eintmodPOTSWITCH)
475 /* We add up the switch cost separately */
476 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_PSW + (calcEnergy ? 1 : 0),
477 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
479 if (ic.vdwtype == evdwPME)
481 /* We add up the LJ Ewald cost separately */
482 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_EWALD + (calcEnergy ? 1 : 0),
483 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
488 nonbonded_verlet_t::dispatchNonbondedKernel(Nbnxm::InteractionLocality iLocality,
489 const interaction_const_t &ic,
493 gmx_enerdata_t *enerd,
496 const PairlistSet &pairlistSet = pairlistSets().pairlistSet(iLocality);
498 switch (kernelSetup().kernelType)
500 case Nbnxm::KernelType::Cpu4x4_PlainC:
501 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
502 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
503 nbnxn_kernel_cpu(pairlistSet,
511 enerd->grpp.ener[egCOULSR],
513 enerd->grpp.ener[egBHAMSR] :
514 enerd->grpp.ener[egLJSR]);
517 case Nbnxm::KernelType::Gpu8x8x8:
518 Nbnxm::gpu_launch_kernel(gpu_nbv, forceFlags, iLocality);
521 case Nbnxm::KernelType::Cpu8x8x8_PlainC:
522 nbnxn_kernel_gpu_ref(pairlistSet.gpuList(),
529 enerd->grpp.ener[egCOULSR],
531 enerd->grpp.ener[egBHAMSR] :
532 enerd->grpp.ener[egLJSR]);
536 GMX_RELEASE_ASSERT(false, "Invalid nonbonded kernel type passed!");
540 accountFlops(nrnb, pairlistSet, *this, ic, forceFlags);
544 nonbonded_verlet_t::dispatchFreeEnergyKernel(Nbnxm::InteractionLocality iLocality,
548 const t_mdatoms &mdatoms,
551 gmx_enerdata_t *enerd,
552 const int forceFlags,
555 const gmx::ArrayRef<t_nblist const * const > nbl_fep = pairlistSets().pairlistSet(iLocality).fepLists();
557 /* When the first list is empty, all are empty and there is nothing to do */
558 if (!pairlistSets().params().haveFep || nbl_fep[0]->nrj == 0)
564 /* Add short-range interactions */
565 donb_flags |= GMX_NONBONDED_DO_SR;
567 /* Currently all group scheme kernels always calculate (shift-)forces */
568 if (forceFlags & GMX_FORCE_FORCES)
570 donb_flags |= GMX_NONBONDED_DO_FORCE;
572 if (forceFlags & GMX_FORCE_VIRIAL)
574 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
576 if (forceFlags & GMX_FORCE_ENERGY)
578 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
581 nb_kernel_data_t kernel_data;
582 real dvdl_nb[efptNR] = { 0 };
583 kernel_data.flags = donb_flags;
584 kernel_data.lambda = lambda;
585 kernel_data.dvdl = dvdl_nb;
587 kernel_data.energygrp_elec = enerd->grpp.ener[egCOULSR];
588 kernel_data.energygrp_vdw = enerd->grpp.ener[egLJSR];
590 GMX_ASSERT(gmx_omp_nthreads_get(emntNonbonded) == nbl_fep.ssize(), "Number of lists should be same as number of NB threads");
592 #pragma omp parallel for schedule(static) num_threads(nbl_fep.ssize())
593 for (int th = 0; th < nbl_fep.ssize(); th++)
597 gmx_nb_free_energy_kernel(nbl_fep[th],
598 x, f, fr, &mdatoms, &kernel_data, nrnb);
600 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
603 if (fepvals->sc_alpha != 0)
605 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
606 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
610 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
611 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
614 /* If we do foreign lambda and we have soft-core interactions
615 * we have to recalculate the (non-linear) energies contributions.
617 if (fepvals->n_lambda > 0 && (forceFlags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
620 kernel_data.flags = (donb_flags & ~(GMX_NONBONDED_DO_FORCE | GMX_NONBONDED_DO_SHIFTFORCE)) | GMX_NONBONDED_DO_FOREIGNLAMBDA;
621 kernel_data.lambda = lam_i;
622 kernel_data.energygrp_elec = enerd->foreign_grpp.ener[egCOULSR];
623 kernel_data.energygrp_vdw = enerd->foreign_grpp.ener[egLJSR];
624 /* Note that we add to kernel_data.dvdl, but ignore the result */
626 for (int i = 0; i < enerd->n_lambda; i++)
628 for (int j = 0; j < efptNR; j++)
630 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
632 reset_foreign_enerdata(enerd);
633 #pragma omp parallel for schedule(static) num_threads(nbl_fep.ssize())
634 for (int th = 0; th < nbl_fep.ssize(); th++)
638 gmx_nb_free_energy_kernel(nbl_fep[th],
639 x, f, fr, &mdatoms, &kernel_data, nrnb);
641 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
644 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
645 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];