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 Translation layer to GROMACS data structures for 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 "nblib/gmxsetup.h"
44 #include "gromacs/ewald/ewald_utils.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/nbnxm.h"
54 #include "gromacs/nbnxm/nbnxm_simd.h"
55 #include "gromacs/nbnxm/pairlistset.h"
56 #include "gromacs/nbnxm/pairlistsets.h"
57 #include "gromacs/nbnxm/pairsearch.h"
58 #include "gromacs/pbcutil/pbc.h"
59 #include "gromacs/utility/logger.h"
60 #include "gromacs/utility/listoflists.h"
61 #include "gromacs/utility/smalloc.h"
62 #include "nblib/exception.h"
63 #include "nblib/kerneloptions.h"
64 #include "nblib/particletype.h"
65 #include "nblib/simulationstate.h"
70 //! Helper to translate between the different enumeration values.
71 static Nbnxm::KernelType translateBenchmarkEnum(const SimdKernels& kernel)
73 int kernelInt = static_cast<int>(kernel);
74 return static_cast<Nbnxm::KernelType>(kernelInt);
77 /*! \brief Checks the kernel setup
79 * Returns an error string when the kernel is not available.
81 static void checkKernelSetup(const NBKernelOptions& options)
83 if (options.nbnxmSimd >= SimdKernels::Count || options.nbnxmSimd == SimdKernels::SimdAuto)
85 throw InputException("Need a valid kernel SIMD type");
88 if ((options.nbnxmSimd != SimdKernels::SimdNo && !GMX_SIMD)
89 #ifndef GMX_NBNXN_SIMD_4XN
90 || options.nbnxmSimd == SimdKernels::Simd4XM
92 #ifndef GMX_NBNXN_SIMD_2XNN
93 || options.nbnxmSimd == SimdKernels::Simd2XMM
97 throw InputException("The requested SIMD kernel was not set up at configuration time");
101 NbvSetupUtil::NbvSetupUtil() : gmxForceCalculator_(std::make_unique<GmxForceCalculator>()) {}
103 void NbvSetupUtil::setExecutionContext(const NBKernelOptions& options)
105 // Todo: find a more general way to initialize hardware
106 gmx_omp_nthreads_set(emntPairsearch, options.numOpenMPThreads);
107 gmx_omp_nthreads_set(emntNonbonded, options.numOpenMPThreads);
110 Nbnxm::KernelSetup NbvSetupUtil::getKernelSetup(const NBKernelOptions& options)
112 checkKernelSetup(options);
114 Nbnxm::KernelSetup kernelSetup;
116 // The int enum options.nbnxnSimd is set up to match Nbnxm::KernelType + 1
117 kernelSetup.kernelType = translateBenchmarkEnum(options.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 = options.useTabulatedEwaldCorr
126 ? Nbnxm::EwaldExclusionType::Table
127 : Nbnxm::EwaldExclusionType::Analytical;
133 void NbvSetupUtil::setParticleInfoAllVdv(const size_t numParticles)
136 particleInfoAllVdw_.resize(numParticles);
137 for (size_t particleI = 0; particleI < numParticles; particleI++)
139 SET_CGINFO_HAS_VDW(particleInfoAllVdw_[particleI]);
140 SET_CGINFO_HAS_Q(particleInfoAllVdw_[particleI]);
144 void NbvSetupUtil::setNonBondedParameters(const std::vector<ParticleType>& particleTypes,
145 const NonBondedInteractionMap& nonBondedInteractionMap)
147 /* Todo: Refactor nbnxm to take nonbondedParameters_ directly
149 * initial self-handling of combination rules
150 * size: 2*(numParticleTypes^2)
152 nonbondedParameters_.reserve(2 * particleTypes.size() * particleTypes.size());
154 constexpr real c6factor = 6.0;
155 constexpr real c12factor = 12.0;
157 for (const ParticleType& particleType1 : particleTypes)
159 for (const ParticleType& particleType2 : particleTypes)
161 nonbondedParameters_.push_back(
162 nonBondedInteractionMap.getC6(particleType1.name(), particleType2.name()) * c6factor);
163 nonbondedParameters_.push_back(
164 nonBondedInteractionMap.getC12(particleType1.name(), particleType2.name()) * c12factor);
169 void NbvSetupUtil::setAtomProperties(const std::vector<int>& particleTypeIdOfAllParticles,
170 const std::vector<real>& charges)
172 gmxForceCalculator_->nbv_->setAtomProperties(particleTypeIdOfAllParticles, charges, particleInfoAllVdw_);
175 //! Sets up and returns a Nbnxm object for the given options and system
176 void NbvSetupUtil::setupNbnxmInstance(const size_t numParticleTypes, const NBKernelOptions& options)
178 const auto pinPolicy = (options.useGpu ? gmx::PinningPolicy::PinnedIfSupported
179 : gmx::PinningPolicy::CannotBePinned);
180 const int numThreads = options.numOpenMPThreads;
181 // Note: the options and Nbnxm combination rule enums values should match
182 const int combinationRule = static_cast<int>(options.ljCombinationRule);
184 checkKernelSetup(options); // throws exception is setup is invalid
186 Nbnxm::KernelSetup kernelSetup = getKernelSetup(options);
188 PairlistParams pairlistParams(kernelSetup.kernelType, false, options.pairlistCutoff, false);
189 Nbnxm::GridSet gridSet(
190 PbcType::Xyz, false, nullptr, nullptr, pairlistParams.pairlistType, false, numThreads, pinPolicy);
191 auto pairlistSets = std::make_unique<PairlistSets>(pairlistParams, false, 0);
192 auto pairSearch = std::make_unique<PairSearch>(
193 PbcType::Xyz, false, nullptr, nullptr, pairlistParams.pairlistType, false, numThreads, pinPolicy);
195 auto atomData = std::make_unique<nbnxn_atomdata_t>(pinPolicy);
197 // Needs to be called with the number of unique ParticleTypes
198 nbnxn_atomdata_init(gmx::MDLogger(),
200 kernelSetup.kernelType,
203 nonbondedParameters_,
207 // Put everything together
208 auto nbv = std::make_unique<nonbonded_verlet_t>(
209 std::move(pairlistSets), std::move(pairSearch), std::move(atomData), kernelSetup, nullptr, nullWallcycle);
211 gmxForceCalculator_->nbv_ = std::move(nbv);
214 //! Computes the Ewald splitting coefficient for Coulomb
215 static real ewaldCoeff(const real ewald_rtol, const real pairlistCutoff)
217 return calc_ewaldcoeff_q(pairlistCutoff, ewald_rtol);
220 void NbvSetupUtil::setupStepWorkload(const NBKernelOptions& options)
222 gmxForceCalculator_->stepWork_->computeForces = true;
223 gmxForceCalculator_->stepWork_->computeNonbondedForces = true;
225 if (options.computeVirialAndEnergy)
227 gmxForceCalculator_->stepWork_->computeVirial = true;
228 gmxForceCalculator_->stepWork_->computeEnergy = true;
232 void NbvSetupUtil::setupInteractionConst(const NBKernelOptions& options)
234 gmxForceCalculator_->interactionConst_->vdwtype = VanDerWaalsType::Cut;
235 gmxForceCalculator_->interactionConst_->vdw_modifier = InteractionModifiers::PotShift;
236 gmxForceCalculator_->interactionConst_->rvdw = options.pairlistCutoff;
238 switch (options.coulombType)
240 case CoulombType::Pme:
241 gmxForceCalculator_->interactionConst_->eeltype = CoulombInteractionType::Pme;
243 case CoulombType::Cutoff:
244 gmxForceCalculator_->interactionConst_->eeltype = CoulombInteractionType::Cut;
246 case CoulombType::ReactionField:
247 gmxForceCalculator_->interactionConst_->eeltype = CoulombInteractionType::RF;
249 case CoulombType::Count: throw InputException("Unsupported electrostatic interaction");
251 gmxForceCalculator_->interactionConst_->coulomb_modifier = InteractionModifiers::PotShift;
252 gmxForceCalculator_->interactionConst_->rcoulomb = options.pairlistCutoff;
253 // Note: values correspond to ic->coulomb_modifier = InteractionModifiers::PotShift
254 gmxForceCalculator_->interactionConst_->dispersion_shift.cpot =
255 -1.0 / gmx::power6(gmxForceCalculator_->interactionConst_->rvdw);
256 gmxForceCalculator_->interactionConst_->repulsion_shift.cpot =
257 -1.0 / gmx::power12(gmxForceCalculator_->interactionConst_->rvdw);
259 // These are the initialized values but we leave them here so that later
260 // these can become options.
261 gmxForceCalculator_->interactionConst_->epsilon_r = 1.0;
262 gmxForceCalculator_->interactionConst_->reactionFieldPermitivity = 1.0;
264 /* Set the Coulomb energy conversion factor */
265 if (gmxForceCalculator_->interactionConst_->epsilon_r != 0)
267 gmxForceCalculator_->interactionConst_->epsfac =
268 ONE_4PI_EPS0 / gmxForceCalculator_->interactionConst_->epsilon_r;
272 /* eps = 0 is infinite dieletric: no Coulomb interactions */
273 gmxForceCalculator_->interactionConst_->epsfac = 0;
277 gmxForceCalculator_->interactionConst_->epsilon_r,
278 gmxForceCalculator_->interactionConst_->reactionFieldPermitivity,
279 gmxForceCalculator_->interactionConst_->rcoulomb,
280 &gmxForceCalculator_->interactionConst_->reactionFieldCoefficient,
281 &gmxForceCalculator_->interactionConst_->reactionFieldShift);
283 if (EEL_PME_EWALD(gmxForceCalculator_->interactionConst_->eeltype))
285 // Ewald coefficients, we ignore the potential shift
286 gmxForceCalculator_->interactionConst_->ewaldcoeff_q = ewaldCoeff(1e-5, options.pairlistCutoff);
287 if (gmxForceCalculator_->interactionConst_->ewaldcoeff_q <= 0)
289 throw InputException("Ewald coefficient should be > 0");
291 gmxForceCalculator_->interactionConst_->coulombEwaldTables =
292 std::make_unique<EwaldCorrectionTables>();
293 init_interaction_const_tables(nullptr, gmxForceCalculator_->interactionConst_.get(), 0, 0);
297 void NbvSetupUtil::setupForceRec(const matrix& box)
299 assert((gmxForceCalculator_->forcerec_ && "Forcerec not initialized"));
300 gmxForceCalculator_->forcerec_->nbfp = nonbondedParameters_;
301 snew(gmxForceCalculator_->forcerec_->shift_vec, numShiftVectors);
302 calc_shifts(box, gmxForceCalculator_->forcerec_->shift_vec);
305 void NbvSetupUtil::setParticlesOnGrid(const std::vector<Vec3>& coordinates, const Box& box)
307 gmxForceCalculator_->setParticlesOnGrid(particleInfoAllVdw_, coordinates, box);
310 void NbvSetupUtil::constructPairList(ExclusionLists<int> exclusionLists)
312 gmx::ListOfLists<int> exclusions(std::move(exclusionLists.ListRanges),
313 std::move(exclusionLists.ListElements));
314 gmxForceCalculator_->nbv_->constructPairlist(
315 gmx::InteractionLocality::Local, exclusions, 0, gmxForceCalculator_->nrnb_.get());
319 std::unique_ptr<GmxForceCalculator> GmxSetupDirector::setupGmxForceCalculator(const SimulationState& system,
320 const NBKernelOptions& options)
322 NbvSetupUtil nbvSetupUtil;
323 nbvSetupUtil.setExecutionContext(options);
324 nbvSetupUtil.setNonBondedParameters(system.topology().getParticleTypes(),
325 system.topology().getNonBondedInteractionMap());
326 nbvSetupUtil.setParticleInfoAllVdv(system.topology().numParticles());
328 nbvSetupUtil.setupInteractionConst(options);
329 nbvSetupUtil.setupStepWorkload(options);
330 nbvSetupUtil.setupNbnxmInstance(system.topology().getParticleTypes().size(), options);
331 nbvSetupUtil.setParticlesOnGrid(system.coordinates(), system.box());
332 nbvSetupUtil.constructPairList(system.topology().exclusionLists());
333 nbvSetupUtil.setAtomProperties(system.topology().getParticleTypeIdOfAllParticles(),
334 system.topology().getCharges());
335 nbvSetupUtil.setupForceRec(system.box().legacyMatrix());
337 return nbvSetupUtil.getGmxForceCalculator();