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/gmx_omp_nthreads.h"
46 #include "gromacs/mdtypes/enerdata.h"
47 #include "gromacs/mdtypes/forceoutput.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/mdtypes/simulation_workload.h"
53 #include "gromacs/nbnxm/gpu_data_mgmt.h"
54 #include "gromacs/nbnxm/nbnxm.h"
55 #include "gromacs/simd/simd.h"
56 #include "gromacs/timing/wallcycle.h"
57 #include "gromacs/utility/gmxassert.h"
58 #include "gromacs/utility/real.h"
60 #include "kernel_common.h"
61 #include "nbnxm_simd.h"
62 #include "pairlistset.h"
63 #include "pairlistsets.h"
64 #include "kernels_reference/kernel_gpu_ref.h"
65 #define INCLUDE_KERNELFUNCTION_TABLES
66 #include "kernels_reference/kernel_ref.h"
67 #ifdef GMX_NBNXN_SIMD_2XNN
68 # include "kernels_simd_2xmm/kernels.h"
70 #ifdef GMX_NBNXN_SIMD_4XN
71 # include "kernels_simd_4xm/kernels.h"
73 #undef INCLUDE_FUNCTION_TABLES
75 /*! \brief Clears the energy group output buffers
77 * \param[in,out] out nbnxn kernel output struct
79 static void clearGroupEnergies(nbnxn_atomdata_output_t* out)
81 std::fill(out->Vvdw.begin(), out->Vvdw.end(), 0.0_real);
82 std::fill(out->Vc.begin(), out->Vc.end(), 0.0_real);
83 std::fill(out->VSvdw.begin(), out->VSvdw.end(), 0.0_real);
84 std::fill(out->VSc.begin(), out->VSc.end(), 0.0_real);
87 /*! \brief Reduce the group-pair energy buffers produced by a SIMD kernel
88 * to single terms in the output buffers.
90 * The SIMD kernels produce a large number of energy buffer in SIMD registers
91 * to avoid scattered reads and writes.
93 * \tparam unrollj The unroll size for j-particles in the SIMD kernel
94 * \param[in] numGroups The number of energy groups
95 * \param[in] numGroups_2log Log2 of numGroups, rounded up
96 * \param[in,out] out Struct with energy buffers
99 static void reduceGroupEnergySimdBuffers(int numGroups, int numGroups_2log, 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] stepWork 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
148 * \param[in] wcycle Pointer to cycle counting data structure.
150 static void nbnxn_kernel_cpu(const PairlistSet& pairlistSet,
151 const Nbnxm::KernelSetup& kernelSetup,
152 nbnxn_atomdata_t* nbat,
153 const interaction_const_t& ic,
155 const gmx::StepWorkload& stepWork,
159 gmx_wallcycle* wcycle)
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;
207 default: GMX_RELEASE_ASSERT(false, "Unknown combination rule");
210 case eintmodFORCESWITCH: vdwkt = vdwktLJFORCESWITCH; break;
211 case eintmodPOTSWITCH: vdwkt = vdwktLJPOTSWITCH; break;
212 default: GMX_RELEASE_ASSERT(false, "Unsupported VdW interaction modifier");
215 else if (ic.vdwtype == evdwPME)
217 if (ic.ljpme_comb_rule == eljpmeGEOM)
219 vdwkt = vdwktLJEWALDCOMBGEOM;
223 vdwkt = vdwktLJEWALDCOMBLB;
224 /* At setup we (should have) selected the C reference kernel */
225 GMX_RELEASE_ASSERT(kernelSetup.kernelType == Nbnxm::KernelType::Cpu4x4_PlainC,
226 "Only the C reference nbnxn SIMD kernel supports LJ-PME with LB "
227 "combination rules");
232 GMX_RELEASE_ASSERT(false, "Unsupported VdW interaction type");
235 gmx::ArrayRef<const NbnxnPairlistCpu> pairlists = pairlistSet.cpuLists();
237 int gmx_unused nthreads = gmx_omp_nthreads_get(emntNonbonded);
238 wallcycle_sub_start(wcycle, ewcsNONBONDED_CLEAR);
239 #pragma omp parallel for schedule(static) num_threads(nthreads)
240 for (gmx::index nb = 0; nb < pairlists.ssize(); nb++)
242 // Presently, the kernels do not call C++ code that can throw,
243 // so no need for a try/catch pair in this OpenMP region.
244 nbnxn_atomdata_output_t* out = &nbat->out[nb];
246 if (clearF == enbvClearFYes)
248 clearForceBuffer(nbat, nb);
250 clear_fshift(out->fshift.data());
255 wallcycle_sub_stop(wcycle, ewcsNONBONDED_CLEAR);
256 wallcycle_sub_start(wcycle, ewcsNONBONDED_KERNEL);
259 // TODO: Change to reference
260 const NbnxnPairlistCpu* pairlist = &pairlists[nb];
262 if (!stepWork.computeEnergy)
264 /* Don't calculate energies */
265 switch (kernelSetup.kernelType)
267 case Nbnxm::KernelType::Cpu4x4_PlainC:
268 nbnxn_kernel_noener_ref[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
270 #ifdef GMX_NBNXN_SIMD_2XNN
271 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
272 nbnxm_kernel_noener_simd_2xmm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
275 #ifdef GMX_NBNXN_SIMD_4XN
276 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
277 nbnxm_kernel_noener_simd_4xm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
280 default: GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
283 else if (out->Vvdw.size() == 1)
285 /* A single energy group (pair) */
289 switch (kernelSetup.kernelType)
291 case Nbnxm::KernelType::Cpu4x4_PlainC:
292 nbnxn_kernel_ener_ref[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
294 #ifdef GMX_NBNXN_SIMD_2XNN
295 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
296 nbnxm_kernel_ener_simd_2xmm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
299 #ifdef GMX_NBNXN_SIMD_4XN
300 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
301 nbnxm_kernel_ener_simd_4xm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
304 default: GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
309 /* Calculate energy group contributions */
310 clearGroupEnergies(out);
314 switch (kernelSetup.kernelType)
316 case Nbnxm::KernelType::Cpu4x4_PlainC:
317 unrollj = c_nbnxnCpuIClusterSize;
318 nbnxn_kernel_energrp_ref[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
320 #ifdef GMX_NBNXN_SIMD_2XNN
321 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
322 unrollj = GMX_SIMD_REAL_WIDTH / 2;
323 nbnxm_kernel_energrp_simd_2xmm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
326 #ifdef GMX_NBNXN_SIMD_4XN
327 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
328 unrollj = GMX_SIMD_REAL_WIDTH;
329 nbnxm_kernel_energrp_simd_4xm[coulkt][vdwkt](pairlist, nbat, &ic, shiftVectors, out);
332 default: GMX_RELEASE_ASSERT(false, "Unsupported kernel architecture");
335 if (kernelSetup.kernelType != Nbnxm::KernelType::Cpu4x4_PlainC)
340 reduceGroupEnergySimdBuffers<2>(nbatParams.nenergrp, nbatParams.neg_2log, out);
343 reduceGroupEnergySimdBuffers<4>(nbatParams.nenergrp, nbatParams.neg_2log, out);
346 reduceGroupEnergySimdBuffers<8>(nbatParams.nenergrp, nbatParams.neg_2log, out);
348 default: GMX_RELEASE_ASSERT(false, "Unsupported j-unroll size");
353 wallcycle_sub_stop(wcycle, ewcsNONBONDED_KERNEL);
355 if (stepWork.computeEnergy)
357 reduce_energies_over_lists(nbat, pairlists.ssize(), vVdw, vCoulomb);
361 static void accountFlops(t_nrnb* nrnb,
362 const PairlistSet& pairlistSet,
363 const nonbonded_verlet_t& nbv,
364 const interaction_const_t& ic,
365 const gmx::StepWorkload& stepWork)
367 const bool usingGpuKernels = nbv.useGpu();
369 int enr_nbnxn_kernel_ljc;
370 if (EEL_RF(ic.eeltype) || ic.eeltype == eelCUT)
372 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_RF;
374 else if ((!usingGpuKernels && nbv.kernelSetup().ewaldExclusionType == Nbnxm::EwaldExclusionType::Analytical)
375 || (usingGpuKernels && Nbnxm::gpu_is_kernel_ewald_analytical(nbv.gpu_nbv)))
377 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_EWALD;
381 enr_nbnxn_kernel_ljc = eNR_NBNXN_LJ_TAB;
383 int enr_nbnxn_kernel_lj = eNR_NBNXN_LJ;
384 if (stepWork.computeEnergy)
386 /* In eNR_??? the nbnxn F+E kernels are always the F kernel + 1 */
387 enr_nbnxn_kernel_ljc += 1;
388 enr_nbnxn_kernel_lj += 1;
391 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc, pairlistSet.natpair_ljq_);
392 inc_nrnb(nrnb, enr_nbnxn_kernel_lj, pairlistSet.natpair_lj_);
393 /* The Coulomb-only kernels are offset -eNR_NBNXN_LJ_RF+eNR_NBNXN_RF */
394 inc_nrnb(nrnb, enr_nbnxn_kernel_ljc - eNR_NBNXN_LJ_RF + eNR_NBNXN_RF, pairlistSet.natpair_q_);
396 if (ic.vdw_modifier == eintmodFORCESWITCH)
398 /* We add up the switch cost separately */
399 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_FSW + (stepWork.computeEnergy ? 1 : 0),
400 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
402 if (ic.vdw_modifier == eintmodPOTSWITCH)
404 /* We add up the switch cost separately */
405 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_PSW + (stepWork.computeEnergy ? 1 : 0),
406 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
408 if (ic.vdwtype == evdwPME)
410 /* We add up the LJ Ewald cost separately */
411 inc_nrnb(nrnb, eNR_NBNXN_ADD_LJ_EWALD + (stepWork.computeEnergy ? 1 : 0),
412 pairlistSet.natpair_ljq_ + pairlistSet.natpair_lj_);
416 void nonbonded_verlet_t::dispatchNonbondedKernel(gmx::InteractionLocality iLocality,
417 const interaction_const_t& ic,
418 const gmx::StepWorkload& stepWork,
420 const t_forcerec& fr,
421 gmx_enerdata_t* enerd,
424 const PairlistSet& pairlistSet = pairlistSets().pairlistSet(iLocality);
426 switch (kernelSetup().kernelType)
428 case Nbnxm::KernelType::Cpu4x4_PlainC:
429 case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
430 case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
431 nbnxn_kernel_cpu(pairlistSet, kernelSetup(), nbat.get(), ic, fr.shift_vec, stepWork,
432 clearF, enerd->grpp.ener[egCOULSR].data(),
433 fr.bBHAM ? enerd->grpp.ener[egBHAMSR].data() : enerd->grpp.ener[egLJSR].data(),
437 case Nbnxm::KernelType::Gpu8x8x8:
438 Nbnxm::gpu_launch_kernel(gpu_nbv, stepWork, iLocality);
441 case Nbnxm::KernelType::Cpu8x8x8_PlainC:
442 nbnxn_kernel_gpu_ref(
443 pairlistSet.gpuList(), nbat.get(), &ic, fr.shift_vec, stepWork, clearF,
444 nbat->out[0].f, nbat->out[0].fshift.data(), enerd->grpp.ener[egCOULSR].data(),
445 fr.bBHAM ? enerd->grpp.ener[egBHAMSR].data() : enerd->grpp.ener[egLJSR].data());
448 default: GMX_RELEASE_ASSERT(false, "Invalid nonbonded kernel type passed!");
451 accountFlops(nrnb, pairlistSet, *this, ic, stepWork);
454 void nonbonded_verlet_t::dispatchFreeEnergyKernel(gmx::InteractionLocality iLocality,
455 const t_forcerec* fr,
457 gmx::ForceWithShiftForces* forceWithShiftForces,
458 const t_mdatoms& mdatoms,
461 gmx_enerdata_t* enerd,
462 const gmx::StepWorkload& stepWork,
465 const auto nbl_fep = pairlistSets().pairlistSet(iLocality).fepLists();
467 /* When the first list is empty, all are empty and there is nothing to do */
468 if (!pairlistSets().params().haveFep || nbl_fep[0]->nrj == 0)
474 /* Add short-range interactions */
475 donb_flags |= GMX_NONBONDED_DO_SR;
477 if (stepWork.computeForces)
479 donb_flags |= GMX_NONBONDED_DO_FORCE;
481 if (stepWork.computeVirial)
483 donb_flags |= GMX_NONBONDED_DO_SHIFTFORCE;
485 if (stepWork.computeEnergy)
487 donb_flags |= GMX_NONBONDED_DO_POTENTIAL;
490 nb_kernel_data_t kernel_data;
491 real dvdl_nb[efptNR] = { 0 };
492 kernel_data.flags = donb_flags;
493 kernel_data.lambda = lambda;
494 kernel_data.dvdl = dvdl_nb;
496 kernel_data.energygrp_elec = enerd->grpp.ener[egCOULSR].data();
497 kernel_data.energygrp_vdw = enerd->grpp.ener[egLJSR].data();
499 GMX_ASSERT(gmx_omp_nthreads_get(emntNonbonded) == nbl_fep.ssize(),
500 "Number of lists should be same as number of NB threads");
502 wallcycle_sub_start(wcycle_, ewcsNONBONDED_FEP);
503 #pragma omp parallel for schedule(static) num_threads(nbl_fep.ssize())
504 for (gmx::index th = 0; th < nbl_fep.ssize(); th++)
508 gmx_nb_free_energy_kernel(nbl_fep[th].get(), x, forceWithShiftForces, fr, &mdatoms,
511 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
514 if (fepvals->sc_alpha != 0)
516 enerd->dvdl_nonlin[efptVDW] += dvdl_nb[efptVDW];
517 enerd->dvdl_nonlin[efptCOUL] += dvdl_nb[efptCOUL];
521 enerd->dvdl_lin[efptVDW] += dvdl_nb[efptVDW];
522 enerd->dvdl_lin[efptCOUL] += dvdl_nb[efptCOUL];
525 /* If we do foreign lambda and we have soft-core interactions
526 * we have to recalculate the (non-linear) energies contributions.
528 if (fepvals->n_lambda > 0 && stepWork.computeDhdl && fepvals->sc_alpha != 0)
531 kernel_data.flags = (donb_flags & ~(GMX_NONBONDED_DO_FORCE | GMX_NONBONDED_DO_SHIFTFORCE))
532 | GMX_NONBONDED_DO_FOREIGNLAMBDA;
533 kernel_data.lambda = lam_i;
534 kernel_data.energygrp_elec = enerd->foreign_grpp.ener[egCOULSR].data();
535 kernel_data.energygrp_vdw = enerd->foreign_grpp.ener[egLJSR].data();
536 /* Note that we add to kernel_data.dvdl, but ignore the result */
538 for (size_t i = 0; i < enerd->enerpart_lambda.size(); i++)
540 for (int j = 0; j < efptNR; j++)
542 lam_i[j] = (i == 0 ? lambda[j] : fepvals->all_lambda[j][i - 1]);
544 reset_foreign_enerdata(enerd);
545 #pragma omp parallel for schedule(static) num_threads(nbl_fep.ssize())
546 for (gmx::index th = 0; th < nbl_fep.ssize(); th++)
550 gmx_nb_free_energy_kernel(nbl_fep[th].get(), x, forceWithShiftForces, fr,
551 &mdatoms, &kernel_data, nrnb);
553 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
556 sum_epot(&(enerd->foreign_grpp), enerd->foreign_term);
557 enerd->enerpart_lambda[i] += enerd->foreign_term[F_EPOT];
560 wallcycle_sub_stop(wcycle_, ewcsNONBONDED_FEP);