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, 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/inputrec.h"
50 #include "gromacs/mdtypes/interaction_const.h"
51 #include "gromacs/mdtypes/md_enums.h"
52 #include "gromacs/mdtypes/mdatom.h"
53 #include "gromacs/mdtypes/simulation_workload.h"
54 #include "gromacs/nbnxm/gpu_data_mgmt.h"
55 #include "gromacs/nbnxm/nbnxm.h"
56 #include "gromacs/simd/simd.h"
57 #include "gromacs/timing/wallcycle.h"
58 #include "gromacs/utility/gmxassert.h"
59 #include "gromacs/utility/real.h"
61 #include "kernel_common.h"
62 #include "nbnxm_gpu.h"
63 #include "nbnxm_simd.h"
64 #include "pairlistset.h"
65 #include "pairlistsets.h"
66 #include "kernels_reference/kernel_gpu_ref.h"
67 #define INCLUDE_KERNELFUNCTION_TABLES
68 #include "kernels_reference/kernel_ref.h"
69 #ifdef GMX_NBNXN_SIMD_2XNN
70 # include "kernels_simd_2xmm/kernels.h"
72 #ifdef GMX_NBNXN_SIMD_4XN
73 # include "kernels_simd_4xm/kernels.h"
75 #undef INCLUDE_FUNCTION_TABLES
77 /*! \brief Clears the energy group output buffers
79 * \param[in,out] out nbnxn kernel output struct
81 static void clearGroupEnergies(nbnxn_atomdata_output_t* out)
83 std::fill(out->Vvdw.begin(), out->Vvdw.end(), 0.0_real);
84 std::fill(out->Vc.begin(), out->Vc.end(), 0.0_real);
85 std::fill(out->VSvdw.begin(), out->VSvdw.end(), 0.0_real);
86 std::fill(out->VSc.begin(), out->VSc.end(), 0.0_real);
89 /*! \brief Reduce the group-pair energy buffers produced by a SIMD kernel
90 * to single terms in the output buffers.
92 * The SIMD kernels produce a large number of energy buffer in SIMD registers
93 * to avoid scattered reads and writes.
95 * \tparam unrollj The unroll size for j-particles in the SIMD kernel
96 * \param[in] numGroups The number of energy groups
97 * \param[in] numGroups_2log Log2 of numGroups, rounded up
98 * \param[in,out] out Struct with energy buffers
100 template<int unrollj>
101 static void reduceGroupEnergySimdBuffers(int numGroups, int numGroups_2log, nbnxn_atomdata_output_t* out)
103 const int unrollj_half = unrollj / 2;
104 /* Energies are stored in SIMD registers with size 2^numGroups_2log */
105 const int numGroupsStorage = (1 << numGroups_2log);
107 const real* gmx_restrict vVdwSimd = out->VSvdw.data();
108 const real* gmx_restrict vCoulombSimd = out->VSc.data();
109 real* gmx_restrict vVdw = out->Vvdw.data();
110 real* gmx_restrict vCoulomb = out->Vc.data();
112 /* The size of the SIMD energy group buffer array is:
113 * numGroups*numGroups*numGroupsStorage*unrollj_half*simd_width
115 for (int i = 0; i < numGroups; i++)
117 for (int j1 = 0; j1 < numGroups; j1++)
119 for (int j0 = 0; j0 < numGroups; j0++)
121 int c = ((i * numGroups + j1) * numGroupsStorage + j0) * unrollj_half * unrollj;
122 for (int s = 0; s < unrollj_half; s++)
124 vVdw[i * numGroups + j0] += vVdwSimd[c + 0];
125 vVdw[i * numGroups + j1] += vVdwSimd[c + 1];
126 vCoulomb[i * numGroups + j0] += vCoulombSimd[c + 0];
127 vCoulomb[i * numGroups + j1] += vCoulombSimd[c + 1];
135 /*! \brief Dispatches the non-bonded N versus M atom cluster CPU kernels.
137 * OpenMP parallelization is performed within this function.
138 * Energy reduction, but not force and shift force reduction, is performed
139 * within this function.
141 * \param[in] pairlistSet Pairlists with local or non-local interactions to compute
142 * \param[in] kernelSetup The non-bonded kernel setup
143 * \param[in,out] nbat The atomdata for the interactions
144 * \param[in] ic Non-bonded interaction constants
145 * \param[in] shiftVectors The PBC shift vectors
146 * \param[in] stepWork Flags that tell what to compute
147 * \param[in] clearF Enum that tells if to clear the force output buffer
148 * \param[out] vCoulomb Output buffer for Coulomb energies
149 * \param[out] vVdw Output buffer for Van der Waals energies
150 * \param[in] wcycle Pointer to cycle counting data structure.
152 static void nbnxn_kernel_cpu(const PairlistSet& pairlistSet,
153 const Nbnxm::KernelSetup& kernelSetup,
154 nbnxn_atomdata_t* nbat,
155 const interaction_const_t& ic,
157 const gmx::StepWorkload& stepWork,
161 gmx_wallcycle* wcycle)
165 if (EEL_RF(ic.eeltype) || ic.eeltype == eelCUT)
171 if (kernelSetup.ewaldExclusionType == Nbnxm::EwaldExclusionType::Table)
173 if (ic.rcoulomb == ic.rvdw)
179 coulkt = coulktTAB_TWIN;
184 if (ic.rcoulomb == ic.rvdw)
186 coulkt = coulktEWALD;
190 coulkt = coulktEWALD_TWIN;
195 const nbnxn_atomdata_t::Params& nbatParams = nbat->params();
198 if (ic.vdwtype == evdwCUT)
200 switch (ic.vdw_modifier)
203 case eintmodPOTSHIFT:
204 switch (nbatParams.comb_rule)
206 case ljcrGEOM: vdwkt = vdwktLJCUT_COMBGEOM; break;
207 case ljcrLB: vdwkt = vdwktLJCUT_COMBLB; break;
208 case ljcrNONE: vdwkt = vdwktLJCUT_COMBNONE; break;
209 default: GMX_RELEASE_ASSERT(false, "Unknown combination rule");
212 case eintmodFORCESWITCH: vdwkt = vdwktLJFORCESWITCH; break;
213 case eintmodPOTSWITCH: vdwkt = vdwktLJPOTSWITCH; break;
214 default: GMX_RELEASE_ASSERT(false, "Unsupported VdW interaction modifier");
217 else if (ic.vdwtype == evdwPME)
219 if (ic.ljpme_comb_rule == eljpmeGEOM)
221 vdwkt = vdwktLJEWALDCOMBGEOM;
225 vdwkt = vdwktLJEWALDCOMBLB;
226 /* At setup we (should have) selected the C reference kernel */
227 GMX_RELEASE_ASSERT(kernelSetup.kernelType == Nbnxm::KernelType::Cpu4x4_PlainC,
228 "Only the C reference nbnxn SIMD kernel supports LJ-PME with LB "
229 "combination rules");
234 GMX_RELEASE_ASSERT(false, "Unsupported VdW interaction type");
237 gmx::ArrayRef<const NbnxnPairlistCpu> pairlists = pairlistSet.cpuLists();
239 int gmx_unused nthreads = gmx_omp_nthreads_get(emntNonbonded);
240 wallcycle_sub_start(wcycle, ewcsNONBONDED_CLEAR);
241 #pragma omp parallel for schedule(static) num_threads(nthreads)
242 for (gmx::index 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());
257 wallcycle_sub_stop(wcycle, ewcsNONBONDED_CLEAR);
258 wallcycle_sub_start(wcycle, ewcsNONBONDED_KERNEL);
261 // TODO: Change to reference
262 const NbnxnPairlistCpu* pairlist = &pairlists[nb];
264 if (!stepWork.computeEnergy)
266 /* Don't calculate energies */
267 switch (kernelSetup.kernelType)
269 case Nbnxm::KernelType::Cpu4x4_PlainC:
270 nbnxn_kernel_noener_ref[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
272 #ifdef GMX_NBNXN_SIMD_2XNN
273 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
274 nbnxm_kernel_noener_simd_2xmm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
277 #ifdef GMX_NBNXN_SIMD_4XN
278 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
279 nbnxm_kernel_noener_simd_4xm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
282 default: GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
285 else if (out->Vvdw.size() == 1)
287 /* A single energy group (pair) */
291 switch (kernelSetup.kernelType)
293 case Nbnxm::KernelType::Cpu4x4_PlainC:
294 nbnxn_kernel_ener_ref[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
296 #ifdef GMX_NBNXN_SIMD_2XNN
297 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
298 nbnxm_kernel_ener_simd_2xmm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
301 #ifdef GMX_NBNXN_SIMD_4XN
302 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
303 nbnxm_kernel_ener_simd_4xm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
306 default: GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
311 /* Calculate energy group contributions */
312 clearGroupEnergies(out);
316 switch (kernelSetup.kernelType)
318 case Nbnxm::KernelType::Cpu4x4_PlainC:
319 unrollj = c_nbnxnCpuIClusterSize;
320 nbnxn_kernel_energrp_ref[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
322 #ifdef GMX_NBNXN_SIMD_2XNN
323 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
324 unrollj = GMX_SIMD_REAL_WIDTH / 2;
325 nbnxm_kernel_energrp_simd_2xmm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
328 #ifdef GMX_NBNXN_SIMD_4XN
329 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
330 unrollj = GMX_SIMD_REAL_WIDTH;
331 nbnxm_kernel_energrp_simd_4xm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
334 default: GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
337 if (kernelSetup.kernelType != Nbnxm::KernelType::Cpu4x4_PlainC)
342 reduceGroupEnergySimdBuffers<2>(nbatParams.nenergrp, nbatParams.neg_2log, out);
345 reduceGroupEnergySimdBuffers<4>(nbatParams.nenergrp, nbatParams.neg_2log, out);
348 reduceGroupEnergySimdBuffers<8>(nbatParams.nenergrp, nbatParams.neg_2log, out);
350 default: GMX_RELEASE_ASSERT(false, "Unsupported j-unroll size");
355 wallcycle_sub_stop(wcycle, ewcsNONBONDED_KERNEL);
357 if (stepWork.computeEnergy)
359 reduce_energies_over_lists(nbat, pairlists.ssize(), vVdw, vCoulomb);
363 static void accountFlops(t_nrnb* nrnb,
364 const PairlistSet& pairlistSet,
365 const nonbonded_verlet_t& nbv,
366 const interaction_const_t& ic,
367 const gmx::StepWorkload& stepWork)
369 const bool usingGpuKernels = nbv.useGpu();
371 int enr_nbnxn_kernel_ljc;
372 if (EEL_RF(ic.eeltype) || ic.eeltype == eelCUT)
374 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
376 else if ((!usingGpuKernels && nbv.kernelSetup().ewaldExclusionType == Nbnxm::EwaldExclusionType::Analytical)
377 || (usingGpuKernels && Nbnxm::gpu_is_kernel_ewald_analytical(nbv.gpu_nbv)))
379 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
383 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
385 int enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
386 if (stepWork.computeEnergy)
388 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
389 enr_nbnxn_kernel_ljc += 1;
390 enr_nbnxn_kernel_lj += 1;
393 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc, pairlistSet.natpair_ljq_);
394 inc_nrnb(nrnb, enr_nbnxn_kernel_lj, pairlistSet.natpair_lj_);
395 /* The Coulomb-only kernels are offset -eNR_NBNXN_LJ_RF+eNR_NBNXN_RF */
396 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc - eNR_NBNXN_LJ_RF + eNR_NBNXN_RF, pairlistSet.natpair_q_);
398 if (ic.vdw_modifier == eintmodFORCESWITCH)
400 /* We add up the switch cost separately */
401 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_FSW + (stepWork.computeEnergy ? 1 : 0),
402 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
404 if (ic.vdw_modifier == eintmodPOTSWITCH)
406 /* We add up the switch cost separately */
407 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_PSW + (stepWork.computeEnergy ? 1 : 0),
408 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
410 if (ic.vdwtype == evdwPME)
412 /* We add up the LJ Ewald cost separately */
413 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_EWALD + (stepWork.computeEnergy ? 1 : 0),
414 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
418 void nonbonded_verlet_t::dispatchNonbondedKernel(gmx::InteractionLocality iLocality,
419 const interaction_const_t& ic,
420 const gmx::StepWorkload& stepWork,
422 const t_forcerec& fr,
423 gmx_enerdata_t* enerd,
426 const PairlistSet& pairlistSet = pairlistSets().pairlistSet(iLocality);
428 switch (kernelSetup().kernelType)
430 case Nbnxm::KernelType::Cpu4x4_PlainC:
431 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
432 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
433 nbnxn_kernel_cpu(pairlistSet, kernelSetup(), nbat.get(), ic, fr.shift_vec, stepWork,
434 clearF, enerd->grpp.ener[egCOULSR].data(),
435 fr.bBHAM ? enerd->grpp.ener[egBHAMSR].data() : enerd->grpp.ener[egLJSR].data(),
439 case Nbnxm::KernelType::Gpu8x8x8:
440 Nbnxm::gpu_launch_kernel(gpu_nbv, stepWork, iLocality);
443 case Nbnxm::KernelType::Cpu8x8x8_PlainC:
444 nbnxn_kernel_gpu_ref(
445 pairlistSet.gpuList(), nbat.get(), &ic, fr.shift_vec, stepWork, clearF,
446 nbat->out[0].f, nbat->out[0].fshift.data(), enerd->grpp.ener[egCOULSR].data(),
447 fr.bBHAM ? enerd->grpp.ener[egBHAMSR].data() : enerd->grpp.ener[egLJSR].data());
450 default: GMX_RELEASE_ASSERT(false, "Invalid nonbonded kernel type passed!");
453 accountFlops(nrnb, pairlistSet, *this, ic, stepWork);
456 void nonbonded_verlet_t::dispatchFreeEnergyKernel(gmx::InteractionLocality iLocality,
457 const t_forcerec* fr,
459 gmx::ForceWithShiftForces* forceWithShiftForces,
460 const t_mdatoms& mdatoms,
463 gmx_enerdata_t* enerd,
464 const gmx::StepWorkload& stepWork,
467 const auto nbl_fep = pairlistSets().pairlistSet(iLocality).fepLists();
469 /* When the first list is empty, all are empty and there is nothing to do */
470 if (!pairlistSets().params().haveFep || nbl_fep[0]->nrj == 0)
476 /* Add short-range interactions */
477 donb_flags |= GMX_NONBONDED_DO_SR;
479 if (stepWork.computeForces)
481 donb_flags |= GMX_NONBONDED_DO_FORCE;
483 if (stepWork.computeVirial)
485 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
487 if (stepWork.computeEnergy)
489 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
492 nb_kernel_data_t kernel_data;
493 real dvdl_nb[efptNR] = { 0 };
494 kernel_data.flags = donb_flags;
495 kernel_data.lambda = lambda;
496 kernel_data.dvdl = dvdl_nb;
498 kernel_data.energygrp_elec = enerd->grpp.ener[egCOULSR].data();
499 kernel_data.energygrp_vdw = enerd->grpp.ener[egLJSR].data();
501 GMX_ASSERT(gmx_omp_nthreads_get(emntNonbonded) == nbl_fep.ssize(),
502 "Number of lists should be same as number of NB threads");
504 wallcycle_sub_start(wcycle_, ewcsNONBONDED_FEP);
505 #pragma omp parallel for schedule(static) num_threads(nbl_fep.ssize())
506 for (gmx::index th = 0; th < nbl_fep.ssize(); th++)
510 gmx_nb_free_energy_kernel(nbl_fep[th].get(), x, forceWithShiftForces, fr, &mdatoms,
513 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
516 if (fepvals->sc_alpha != 0)
518 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
519 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
523 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
524 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
527 /* If we do foreign lambda and we have soft-core interactions
528 * we have to recalculate the (non-linear) energies contributions.
530 if (fepvals->n_lambda > 0 && stepWork.computeDhdl && fepvals->sc_alpha != 0)
533 kernel_data.flags = (donb_flags & ~(GMX_NONBONDED_DO_FORCE | GMX_NONBONDED_DO_SHIFTFORCE))
534 | GMX_NONBONDED_DO_FOREIGNLAMBDA;
535 kernel_data.lambda = lam_i;
536 kernel_data.energygrp_elec = enerd->foreign_grpp.ener[egCOULSR].data();
537 kernel_data.energygrp_vdw = enerd->foreign_grpp.ener[egLJSR].data();
538 /* Note that we add to kernel_data.dvdl, but ignore the result */
540 for (size_t i = 0; i < enerd->enerpart_lambda.size(); i++)
542 for (int j = 0; j < efptNR; j++)
544 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i - 1]);
546 reset_foreign_enerdata(enerd);
547 #pragma omp parallel for schedule(static) num_threads(nbl_fep.ssize())
548 for (gmx::index th = 0; th < nbl_fep.ssize(); th++)
552 gmx_nb_free_energy_kernel(nbl_fep[th].get(), x, forceWithShiftForces, fr,
553 &mdatoms, &kernel_data, nrnb);
555 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
558 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
559 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
562 wallcycle_sub_stop(wcycle_, ewcsNONBONDED_FEP);