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/pairlist.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 #define INCLUDE_KERNELFUNCTION_TABLES
62 #include "gromacs/nbnxm/kernels_reference/kernel_ref.h"
63 #ifdef GMX_NBNXN_SIMD_2XNN
64 #include "gromacs/nbnxm/kernels_simd_2xmm/kernels.h"
66 #ifdef GMX_NBNXN_SIMD_4XN
67 #include "gromacs/nbnxm/kernels_simd_4xm/kernels.h"
69 #undef INCLUDE_FUNCTION_TABLES
71 /*! \brief Clears the energy group output buffers
73 * \param[in,out] out nbnxn kernel output struct
75 static void clearGroupEnergies(nbnxn_atomdata_output_t *out)
77 std::fill(out->Vvdw.begin(), out->Vvdw.end(), 0.0_real);
78 std::fill(out->Vc.begin(), out->Vc.end(), 0.0_real);
79 std::fill(out->VSvdw.begin(), out->VSvdw.end(), 0.0_real);
80 std::fill(out->VSc.begin(), out->VSc.end(), 0.0_real);
83 /*! \brief Reduce the group-pair energy buffers produced by a SIMD kernel
84 * to single terms in the output buffers.
86 * The SIMD kernels produce a large number of energy buffer in SIMD registers
87 * to avoid scattered reads and writes.
89 * \tparam unrollj The unroll size for j-particles in the SIMD kernel
90 * \param[in] numGroups The number of energy groups
91 * \param[in] numGroups_2log Log2 of numGroups, rounded up
92 * \param[in,out] out Struct with energy buffers
94 template <int unrollj> static void
95 reduceGroupEnergySimdBuffers(int numGroups,
97 nbnxn_atomdata_output_t *out)
99 const int unrollj_half = unrollj/2;
100 /* Energies are stored in SIMD registers with size 2^numGroups_2log */
101 const int numGroupsStorage = (1 << numGroups_2log);
103 const real * gmx_restrict vVdwSimd = out->VSvdw.data();
104 const real * gmx_restrict vCoulombSimd = out->VSc.data();
105 real * gmx_restrict vVdw = out->Vvdw.data();
106 real * gmx_restrict vCoulomb = out->Vc.data();
108 /* The size of the SIMD energy group buffer array is:
109 * numGroups*numGroups*numGroupsStorage*unrollj_half*simd_width
111 for (int i = 0; i < numGroups; i++)
113 for (int j1 = 0; j1 < numGroups; j1++)
115 for (int j0 = 0; j0 < numGroups; j0++)
117 int c = ((i*numGroups + j1)*numGroupsStorage + j0)*unrollj_half*unrollj;
118 for (int s = 0; s < unrollj_half; s++)
120 vVdw [i*numGroups + j0] += vVdwSimd [c + 0];
121 vVdw [i*numGroups + j1] += vVdwSimd [c + 1];
122 vCoulomb[i*numGroups + j0] += vCoulombSimd[c + 0];
123 vCoulomb[i*numGroups + j1] += vCoulombSimd[c + 1];
131 /*! \brief Dispatches the non-bonded N versus M atom cluster CPU kernels.
133 * OpenMP parallelization is performed within this function.
134 * Energy reduction, but not force and shift force reduction, is performed
135 * within this function.
137 * \param[in] pairlistSet Pairlists with local or non-local interactions to compute
138 * \param[in] kernelSetup The non-bonded kernel setup
139 * \param[in,out] nbat The atomdata for the interactions
140 * \param[in] ic Non-bonded interaction constants
141 * \param[in] shiftVectors The PBC shift vectors
142 * \param[in] forceFlags Flags that tell what to compute
143 * \param[in] clearF Enum that tells if to clear the force output buffer
144 * \param[out] fshift Shift force output buffer
145 * \param[out] vCoulomb Output buffer for Coulomb energies
146 * \param[out] vVdw Output buffer for Van der Waals energies
149 nbnxn_kernel_cpu(const PairlistSet &pairlistSet,
150 const Nbnxm::KernelSetup &kernelSetup,
151 nbnxn_atomdata_t *nbat,
152 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 clear_f(nbat, nb, out->f.data());
254 if ((forceFlags & GMX_FORCE_VIRIAL) && pairlists.ssize() == 1)
260 fshift_p = out->fshift.data();
262 if (clearF == enbvClearFYes)
264 clear_fshift(fshift_p);
268 // TODO: Change to reference
269 const NbnxnPairlistCpu *pairlist = &pairlists[nb];
271 if (!(forceFlags & GMX_FORCE_ENERGY))
273 /* Don't calculate energies */
274 switch (kernelSetup.kernelType)
276 case Nbnxm::KernelType::Cpu4x4_PlainC:
277 nbnxn_kernel_noener_ref[coulkt][vdwkt](pairlist, nbat,
283 #ifdef GMX_NBNXN_SIMD_2XNN
284 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
285 nbnxm_kernel_noener_simd_2xmm[coulkt][vdwkt](pairlist, nbat,
292 #ifdef GMX_NBNXN_SIMD_4XN
293 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
294 nbnxm_kernel_noener_simd_4xm[coulkt][vdwkt](pairlist, nbat,
302 GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
305 else if (out->Vvdw.size() == 1)
307 /* A single energy group (pair) */
311 switch (kernelSetup.kernelType)
313 case Nbnxm::KernelType::Cpu4x4_PlainC:
314 nbnxn_kernel_ener_ref[coulkt][vdwkt](pairlist, nbat,
322 #ifdef GMX_NBNXN_SIMD_2XNN
323 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
324 nbnxm_kernel_ener_simd_2xmm[coulkt][vdwkt](pairlist, nbat,
333 #ifdef GMX_NBNXN_SIMD_4XN
334 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
335 nbnxm_kernel_ener_simd_4xm[coulkt][vdwkt](pairlist, nbat,
345 GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
350 /* Calculate energy group contributions */
351 clearGroupEnergies(out);
355 switch (kernelSetup.kernelType)
357 case Nbnxm::KernelType::Cpu4x4_PlainC:
358 unrollj = c_nbnxnCpuIClusterSize;
359 nbnxn_kernel_energrp_ref[coulkt][vdwkt](pairlist, nbat,
367 #ifdef GMX_NBNXN_SIMD_2XNN
368 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
369 unrollj = GMX_SIMD_REAL_WIDTH/2;
370 nbnxm_kernel_energrp_simd_2xmm[coulkt][vdwkt](pairlist, nbat,
379 #ifdef GMX_NBNXN_SIMD_4XN
380 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
381 unrollj = GMX_SIMD_REAL_WIDTH;
382 nbnxm_kernel_energrp_simd_4xm[coulkt][vdwkt](pairlist, nbat,
392 GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
395 if (kernelSetup.kernelType != Nbnxm::KernelType::Cpu4x4_PlainC)
400 reduceGroupEnergySimdBuffers<2>(nbatParams.nenergrp,
405 reduceGroupEnergySimdBuffers<4>(nbatParams.nenergrp,
410 reduceGroupEnergySimdBuffers<8>(nbatParams.nenergrp,
415 GMX_RELEASE_ASSERT(false, "Unsupported j-unroll size");
421 if (forceFlags & GMX_FORCE_ENERGY)
423 reduce_energies_over_lists(nbat, pairlists.ssize(), vVdw, vCoulomb);
427 static void accountFlops(t_nrnb *nrnb,
428 const PairlistSet &pairlistSet,
429 const nonbonded_verlet_t &nbv,
430 const interaction_const_t &ic,
431 const int forceFlags)
433 const bool usingGpuKernels = nbv.useGpu();
435 int enr_nbnxn_kernel_ljc;
436 if (EEL_RF(ic.eeltype) || ic.eeltype == eelCUT)
438 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
440 else if ((!usingGpuKernels && nbv.kernelSetup().ewaldExclusionType == Nbnxm::EwaldExclusionType::Analytical) ||
441 (usingGpuKernels && Nbnxm::gpu_is_kernel_ewald_analytical(nbv.gpu_nbv)))
443 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
447 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
449 int enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
450 if (forceFlags & GMX_FORCE_ENERGY)
452 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
453 enr_nbnxn_kernel_ljc += 1;
454 enr_nbnxn_kernel_lj += 1;
457 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc,
458 pairlistSet.natpair_ljq_);
459 inc_nrnb(nrnb, enr_nbnxn_kernel_lj,
460 pairlistSet.natpair_lj_);
461 /* The Coulomb-only kernels are offset -eNR_NBNXN_LJ_RF+eNR_NBNXN_RF */
462 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc-eNR_NBNXN_LJ_RF+eNR_NBNXN_RF,
463 pairlistSet.natpair_q_);
465 const bool calcEnergy = ((forceFlags & GMX_FORCE_ENERGY) != 0);
466 if (ic.vdw_modifier == eintmodFORCESWITCH)
468 /* We add up the switch cost separately */
469 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_FSW + (calcEnergy ? 1 : 0),
470 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
472 if (ic.vdw_modifier == eintmodPOTSWITCH)
474 /* We add up the switch cost separately */
475 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_PSW + (calcEnergy ? 1 : 0),
476 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
478 if (ic.vdwtype == evdwPME)
480 /* We add up the LJ Ewald cost separately */
481 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_EWALD + (calcEnergy ? 1 : 0),
482 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
487 nonbonded_verlet_t::dispatchNonbondedKernel(Nbnxm::InteractionLocality iLocality,
488 const interaction_const_t &ic,
492 gmx_enerdata_t *enerd,
495 const PairlistSet &pairlistSet = pairlistSets().pairlistSet(iLocality);
497 switch (kernelSetup().kernelType)
499 case Nbnxm::KernelType::Cpu4x4_PlainC:
500 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
501 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
502 nbnxn_kernel_cpu(pairlistSet,
510 enerd->grpp.ener[egCOULSR],
512 enerd->grpp.ener[egBHAMSR] :
513 enerd->grpp.ener[egLJSR]);
516 case Nbnxm::KernelType::Gpu8x8x8:
517 Nbnxm::gpu_launch_kernel(gpu_nbv, forceFlags, iLocality);
520 case Nbnxm::KernelType::Cpu8x8x8_PlainC:
521 nbnxn_kernel_gpu_ref(pairlistSet.gpuList(),
528 enerd->grpp.ener[egCOULSR],
530 enerd->grpp.ener[egBHAMSR] :
531 enerd->grpp.ener[egLJSR]);
535 GMX_RELEASE_ASSERT(false, "Invalid nonbonded kernel type passed!");
539 accountFlops(nrnb, pairlistSet, *this, ic, forceFlags);
543 nonbonded_verlet_t::dispatchFreeEnergyKernel(Nbnxm::InteractionLocality iLocality,
547 const t_mdatoms &mdatoms,
550 gmx_enerdata_t *enerd,
551 const int forceFlags,
554 const gmx::ArrayRef<t_nblist const * const > nbl_fep = pairlistSets().pairlistSet(iLocality).fepLists();
556 /* When the first list is empty, all are empty and there is nothing to do */
557 if (!pairlistSets().params().haveFep || nbl_fep[0]->nrj == 0)
563 /* Add short-range interactions */
564 donb_flags |= GMX_NONBONDED_DO_SR;
566 /* Currently all group scheme kernels always calculate (shift-)forces */
567 if (forceFlags & GMX_FORCE_FORCES)
569 donb_flags |= GMX_NONBONDED_DO_FORCE;
571 if (forceFlags & GMX_FORCE_VIRIAL)
573 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
575 if (forceFlags & GMX_FORCE_ENERGY)
577 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
580 nb_kernel_data_t kernel_data;
581 real dvdl_nb[efptNR] = { 0 };
582 kernel_data.flags = donb_flags;
583 kernel_data.lambda = lambda;
584 kernel_data.dvdl = dvdl_nb;
586 kernel_data.energygrp_elec = enerd->grpp.ener[egCOULSR];
587 kernel_data.energygrp_vdw = enerd->grpp.ener[egLJSR];
589 GMX_ASSERT(gmx_omp_nthreads_get(emntNonbonded) == nbl_fep.ssize(), "Number of lists should be same as number of NB threads");
591 #pragma omp parallel for schedule(static) num_threads(nbl_fep.ssize())
592 for (int th = 0; th < nbl_fep.ssize(); th++)
596 gmx_nb_free_energy_kernel(nbl_fep[th],
597 x, f, fr, &mdatoms, &kernel_data, nrnb);
599 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
602 if (fepvals->sc_alpha != 0)
604 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
605 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
609 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
610 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
613 /* If we do foreign lambda and we have soft-core interactions
614 * we have to recalculate the (non-linear) energies contributions.
616 if (fepvals->n_lambda > 0 && (forceFlags & GMX_FORCE_DHDL) && fepvals->sc_alpha != 0)
619 kernel_data.flags = (donb_flags & ~(GMX_NONBONDED_DO_FORCE | GMX_NONBONDED_DO_SHIFTFORCE)) | GMX_NONBONDED_DO_FOREIGNLAMBDA;
620 kernel_data.lambda = lam_i;
621 kernel_data.energygrp_elec = enerd->foreign_grpp.ener[egCOULSR];
622 kernel_data.energygrp_vdw = enerd->foreign_grpp.ener[egLJSR];
623 /* Note that we add to kernel_data.dvdl, but ignore the result */
625 for (int i = 0; i < enerd->n_lambda; i++)
627 for (int j = 0; j < efptNR; j++)
629 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i-1]);
631 reset_foreign_enerdata(enerd);
632 #pragma omp parallel for schedule(static) num_threads(nbl_fep.ssize())
633 for (int th = 0; th < nbl_fep.ssize(); th++)
637 gmx_nb_free_energy_kernel(nbl_fep[th],
638 x, f, fr, &mdatoms, &kernel_data, nrnb);
640 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
643 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
644 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];