2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2020,2021, 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.
36 * \brief Utilities to setup GROMACS data structures for non-bonded force calculations.
38 * \author Victor Holanda <victor.holanda@cscs.ch>
39 * \author Joe Jordan <ejjordan@kth.se>
40 * \author Prashanth Kanduri <kanduri@cscs.ch>
41 * \author Sebastian Keller <keller@cscs.ch>
43 #include "gromacs/ewald/ewald_utils.h"
44 #include "gromacs/gpu_utils/device_stream_manager.h"
45 #include "gromacs/gmxlib/nrnb.h"
46 #include "gromacs/mdlib/forcerec.h"
47 #include "gromacs/mdlib/gmx_omp_nthreads.h"
48 #include "gromacs/mdlib/rf_util.h"
49 #include "gromacs/mdtypes/forcerec.h"
50 #include "gromacs/mdtypes/interaction_const.h"
51 #include "gromacs/mdtypes/simulation_workload.h"
52 #include "gromacs/nbnxm/atomdata.h"
53 #include "gromacs/nbnxm/gpu_data_mgmt.h"
54 #include "gromacs/nbnxm/nbnxm_gpu.h"
55 #include "gromacs/nbnxm/nbnxm.h"
56 #include "gromacs/nbnxm/nbnxm_simd.h"
57 #include "gromacs/nbnxm/pairlistset.h"
58 #include "gromacs/nbnxm/pairlistsets.h"
59 #include "gromacs/nbnxm/pairsearch.h"
60 #include "gromacs/pbcutil/pbc.h"
61 #include "gromacs/utility/logger.h"
62 #include "gromacs/utility/smalloc.h"
63 #include "nblib/exception.h"
64 #include "nblib/kerneloptions.h"
65 #include "nblib/particletype.h"
67 #include "nbnxmsetuphelpers.h"
72 int64_t findNumEnergyGroups(gmx::ArrayRef<int64_t> particleInteractionFlags)
74 auto groupId = [](int code1, int code2) {
75 return (code1 & gmx::sc_atomInfo_EnergyGroupIdMask) < (code2 & gmx::sc_atomInfo_EnergyGroupIdMask);
78 int maxElement = *std::max_element(
79 std::begin(particleInteractionFlags), std::end(particleInteractionFlags), groupId);
80 return ((maxElement + 1) & gmx::sc_atomInfo_EnergyGroupIdMask);
83 Nbnxm::KernelType translateBenchmarkEnum(const SimdKernels& kernel)
85 int kernelInt = static_cast<int>(kernel);
86 return static_cast<Nbnxm::KernelType>(kernelInt);
89 void checkKernelSetupSimd(const SimdKernels nbnxmSimd)
91 if (nbnxmSimd >= SimdKernels::Count || nbnxmSimd == SimdKernels::SimdAuto)
93 throw InputException("Need a valid kernel SIMD type");
96 if ((nbnxmSimd != SimdKernels::SimdNo && !GMX_SIMD)
97 #ifndef GMX_NBNXN_SIMD_4XN
98 || nbnxmSimd == SimdKernels::Simd4XM
100 #ifndef GMX_NBNXN_SIMD_2XNN
101 || nbnxmSimd == SimdKernels::Simd2XMM
105 throw InputException("The requested SIMD kernel was not set up at configuration time");
109 Nbnxm::KernelSetup createKernelSetupCPU(const SimdKernels nbnxmSimd, const bool useTabulatedEwaldCorr)
111 checkKernelSetupSimd(nbnxmSimd);
113 Nbnxm::KernelSetup kernelSetup;
115 // The int enum options.nbnxnSimd is set up to match Nbnxm::KernelType + 1
116 kernelSetup.kernelType = translateBenchmarkEnum(nbnxmSimd);
118 // The plain-C kernel does not support analytical ewald correction
119 if (kernelSetup.kernelType == Nbnxm::KernelType::Cpu4x4_PlainC)
121 kernelSetup.ewaldExclusionType = Nbnxm::EwaldExclusionType::Table;
125 kernelSetup.ewaldExclusionType = useTabulatedEwaldCorr ? Nbnxm::EwaldExclusionType::Table
126 : Nbnxm::EwaldExclusionType::Analytical;
132 Nbnxm::KernelSetup createKernelSetupGPU(const bool useTabulatedEwaldCorr)
134 Nbnxm::KernelSetup kernelSetup;
135 kernelSetup.kernelType = Nbnxm::KernelType::Gpu8x8x8;
136 kernelSetup.ewaldExclusionType = useTabulatedEwaldCorr ? Nbnxm::EwaldExclusionType::Table
137 : Nbnxm::EwaldExclusionType::Analytical;
142 std::vector<int64_t> createParticleInfoAllVdw(const size_t numParticles)
144 std::vector<int64_t> particleInfoAllVdw(numParticles);
145 for (size_t particleI = 0; particleI < numParticles; particleI++)
147 particleInfoAllVdw[particleI] |= gmx::sc_atomInfo_HasVdw;
148 particleInfoAllVdw[particleI] |= gmx::sc_atomInfo_HasCharge;
150 return particleInfoAllVdw;
153 std::vector<real> createNonBondedParameters(const std::vector<ParticleType>& particleTypes,
154 const NonBondedInteractionMap& nonBondedInteractionMap)
156 /* Todo: Refactor nbnxm to take nonbondedParameters_ directly
158 * initial self-handling of combination rules
159 * size: 2*(numParticleTypes^2)
161 std::vector<real> nonbondedParameters;
162 nonbondedParameters.reserve(2 * particleTypes.size() * particleTypes.size());
164 constexpr real c6factor = 6.0;
165 constexpr real c12factor = 12.0;
167 for (const ParticleType& particleType1 : particleTypes)
169 for (const ParticleType& particleType2 : particleTypes)
171 nonbondedParameters.push_back(
172 nonBondedInteractionMap.getC6(particleType1.name(), particleType2.name()) * c6factor);
173 nonbondedParameters.push_back(
174 nonBondedInteractionMap.getC12(particleType1.name(), particleType2.name()) * c12factor);
177 return nonbondedParameters;
180 gmx::StepWorkload createStepWorkload()
182 gmx::StepWorkload stepWorkload;
183 stepWorkload.computeForces = true;
184 stepWorkload.computeNonbondedForces = true;
185 stepWorkload.useGpuFBufferOps = false;
186 stepWorkload.useGpuXBufferOps = false;
191 static gmx::SimulationWorkload createSimulationWorkload()
193 gmx::SimulationWorkload simulationWork;
194 simulationWork.computeNonbonded = true;
195 return simulationWork;
198 gmx::SimulationWorkload createSimulationWorkloadGpu()
200 gmx::SimulationWorkload simulationWork = createSimulationWorkload();
202 simulationWork.useGpuNonbonded = true;
203 simulationWork.useGpuUpdate = false;
205 return simulationWork;
208 std::shared_ptr<gmx::DeviceStreamManager> createDeviceStreamManager(const DeviceInformation& deviceInfo,
209 const gmx::SimulationWorkload& simulationWorkload)
211 return std::make_shared<gmx::DeviceStreamManager>(deviceInfo, false, simulationWorkload, false);
214 real ewaldCoeff(const real ewald_rtol, const real pairlistCutoff)
216 return calc_ewaldcoeff_q(pairlistCutoff, ewald_rtol);
219 interaction_const_t createInteractionConst(const NBKernelOptions& options)
221 interaction_const_t interactionConst;
222 interactionConst.vdwtype = VanDerWaalsType::Cut;
223 interactionConst.vdw_modifier = InteractionModifiers::PotShift;
224 interactionConst.rvdw = options.pairlistCutoff;
226 switch (options.coulombType)
228 case CoulombType::Pme: interactionConst.eeltype = CoulombInteractionType::Pme; break;
229 case CoulombType::Cutoff: interactionConst.eeltype = CoulombInteractionType::Cut; break;
230 case CoulombType::ReactionField:
231 interactionConst.eeltype = CoulombInteractionType::RF;
233 case CoulombType::Count: throw InputException("Unsupported electrostatic interaction");
235 interactionConst.coulomb_modifier = InteractionModifiers::PotShift;
236 interactionConst.rcoulomb = options.pairlistCutoff;
237 // Note: values correspond to ic->coulomb_modifier = eintmodPOTSHIFT
238 interactionConst.dispersion_shift.cpot = -1.0 / gmx::power6(interactionConst.rvdw);
239 interactionConst.repulsion_shift.cpot = -1.0 / gmx::power12(interactionConst.rvdw);
241 // These are the initialized values but we leave them here so that later
242 // these can become options.
243 interactionConst.epsilon_r = 1.0;
244 interactionConst.reactionFieldPermitivity = 1.0;
246 /* Set the Coulomb energy conversion factor */
247 if (interactionConst.epsilon_r != 0)
249 interactionConst.epsfac = ONE_4PI_EPS0 / interactionConst.epsilon_r;
253 /* eps = 0 is infinite dieletric: no Coulomb interactions */
254 interactionConst.epsfac = 0;
258 interactionConst.epsilon_r,
259 interactionConst.reactionFieldPermitivity,
260 interactionConst.rcoulomb,
261 &interactionConst.reactionFieldCoefficient,
262 &interactionConst.reactionFieldShift);
265 if (EEL_PME_EWALD(interactionConst.eeltype))
267 // Ewald coefficients, we ignore the potential shift
268 interactionConst.ewaldcoeff_q = ewaldCoeff(1e-5, options.pairlistCutoff);
269 if (interactionConst.ewaldcoeff_q <= 0)
271 throw InputException("Ewald coefficient should be > 0");
273 interactionConst.coulombEwaldTables = std::make_unique<EwaldCorrectionTables>();
274 init_interaction_const_tables(nullptr, &interactionConst, 0, 0);
276 return interactionConst;
279 std::unique_ptr<nonbonded_verlet_t> createNbnxmCPU(const size_t numParticleTypes,
280 const NBKernelOptions& options,
282 gmx::ArrayRef<const real> nonbondedParameters)
284 const auto pinPolicy = gmx::PinningPolicy::CannotBePinned;
285 const int numThreads = options.numOpenMPThreads;
286 // Note: the options and Nbnxm combination rule enums values should match
287 const int combinationRule = static_cast<int>(options.ljCombinationRule);
289 Nbnxm::KernelSetup kernelSetup =
290 createKernelSetupCPU(options.nbnxmSimd, options.useTabulatedEwaldCorr);
292 PairlistParams pairlistParams(kernelSetup.kernelType, false, options.pairlistCutoff, false);
294 auto pairlistSets = std::make_unique<PairlistSets>(pairlistParams, false, 0);
295 auto pairSearch = std::make_unique<PairSearch>(
296 PbcType::Xyz, false, nullptr, nullptr, pairlistParams.pairlistType, false, numThreads, pinPolicy);
298 // Needs to be called with the number of unique ParticleTypes
299 auto atomData = std::make_unique<nbnxn_atomdata_t>(pinPolicy,
301 kernelSetup.kernelType,
308 // Put everything together
309 auto nbv = std::make_unique<nonbonded_verlet_t>(
310 std::move(pairlistSets), std::move(pairSearch), std::move(atomData), kernelSetup, nullptr, nullptr);
315 std::unique_ptr<nonbonded_verlet_t> createNbnxmGPU(const size_t numParticleTypes,
316 const NBKernelOptions& options,
317 const std::vector<real>& nonbondedParameters,
318 const interaction_const_t& interactionConst,
319 const gmx::DeviceStreamManager& deviceStreamManager)
321 const auto pinPolicy = gmx::PinningPolicy::PinnedIfSupported;
322 const int combinationRule = static_cast<int>(options.ljCombinationRule);
324 Nbnxm::KernelSetup kernelSetup = createKernelSetupGPU(options.useTabulatedEwaldCorr);
326 PairlistParams pairlistParams(kernelSetup.kernelType, false, options.pairlistCutoff, false);
329 // nbnxn_atomdata is always initialized with 1 thread if the GPU is used
330 constexpr int numThreadsInit = 1;
331 // multiple energy groups are not supported on the GPU
332 constexpr int numEnergyGroups = 1;
333 auto atomData = std::make_unique<nbnxn_atomdata_t>(pinPolicy,
335 kernelSetup.kernelType,
342 NbnxmGpu* nbnxmGpu = Nbnxm::gpu_init(
343 deviceStreamManager, &interactionConst, pairlistParams, atomData.get(), false);
345 // minimum iList count for GPU balancing
346 int iListCount = Nbnxm::gpu_min_ci_balanced(nbnxmGpu);
348 auto pairlistSets = std::make_unique<PairlistSets>(pairlistParams, false, iListCount);
349 auto pairSearch = std::make_unique<PairSearch>(
350 PbcType::Xyz, false, nullptr, nullptr, pairlistParams.pairlistType, false, options.numOpenMPThreads, pinPolicy);
352 // Put everything together
353 auto nbv = std::make_unique<nonbonded_verlet_t>(
354 std::move(pairlistSets), std::move(pairSearch), std::move(atomData), kernelSetup, nbnxmGpu, nullptr);
356 // Some paramters must be copied to NbnxmGpu to have a fully constructed nonbonded_verlet_t
357 Nbnxm::gpu_init_atomdata(nbv->gpu_nbv, nbv->nbat.get());
362 void setGmxNonBondedNThreads(int numThreads)
364 gmx_omp_nthreads_set(ModuleMultiThread::Pairsearch, numThreads);
365 gmx_omp_nthreads_set(ModuleMultiThread::Nonbonded, numThreads);
368 void updateForcerec(t_forcerec* forcerec, const matrix& box)
370 assert(forcerec != nullptr && "Forcerec not initialized");
371 forcerec->shift_vec.resize(numShiftVectors);
372 calc_shifts(box, forcerec->shift_vec);