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/atominfo.h"
50 #include "gromacs/mdtypes/forcerec.h"
51 #include "gromacs/mdtypes/interaction_const.h"
52 #include "gromacs/mdtypes/simulation_workload.h"
53 #include "gromacs/nbnxm/atomdata.h"
54 #include "gromacs/nbnxm/nbnxm.h"
55 #include "gromacs/nbnxm/nbnxm_simd.h"
56 #include "gromacs/nbnxm/pairlistset.h"
57 #include "gromacs/nbnxm/pairlistsets.h"
58 #include "gromacs/nbnxm/pairsearch.h"
59 #include "gromacs/pbcutil/pbc.h"
60 #include "gromacs/utility/logger.h"
61 #include "gromacs/utility/listoflists.h"
62 #include "gromacs/utility/smalloc.h"
63 #include "nblib/exception.h"
64 #include "nblib/kerneloptions.h"
65 #include "nblib/particletype.h"
66 #include "nblib/simulationstate.h"
71 //! Helper to translate between the different enumeration values.
72 static Nbnxm::KernelType translateBenchmarkEnum(const SimdKernels& kernel)
74 int kernelInt = static_cast<int>(kernel);
75 return static_cast<Nbnxm::KernelType>(kernelInt);
78 /*! \brief Checks the kernel setup
80 * Returns an error string when the kernel is not available.
82 static void checkKernelSetup(const NBKernelOptions& options)
84 if (options.nbnxmSimd >= SimdKernels::Count || options.nbnxmSimd == SimdKernels::SimdAuto)
86 throw InputException("Need a valid kernel SIMD type");
89 if ((options.nbnxmSimd != SimdKernels::SimdNo && !GMX_SIMD)
90 #ifndef GMX_NBNXN_SIMD_4XN
91 || options.nbnxmSimd == SimdKernels::Simd4XM
93 #ifndef GMX_NBNXN_SIMD_2XNN
94 || options.nbnxmSimd == SimdKernels::Simd2XMM
98 throw InputException("The requested SIMD kernel was not set up at configuration time");
102 NbvSetupUtil::NbvSetupUtil() : gmxForceCalculator_(std::make_unique<GmxForceCalculator>()) {}
104 void NbvSetupUtil::setExecutionContext(const NBKernelOptions& options)
106 // Todo: find a more general way to initialize hardware
107 gmx_omp_nthreads_set(ModuleMultiThread::Pairsearch, options.numOpenMPThreads);
108 gmx_omp_nthreads_set(ModuleMultiThread::Nonbonded, options.numOpenMPThreads);
111 Nbnxm::KernelSetup NbvSetupUtil::getKernelSetup(const NBKernelOptions& options)
113 checkKernelSetup(options);
115 Nbnxm::KernelSetup kernelSetup;
117 // The int enum options.nbnxnSimd is set up to match Nbnxm::KernelType + 1
118 kernelSetup.kernelType = translateBenchmarkEnum(options.nbnxmSimd);
119 // The plain-C kernel does not support analytical ewald correction
120 if (kernelSetup.kernelType == Nbnxm::KernelType::Cpu4x4_PlainC)
122 kernelSetup.ewaldExclusionType = Nbnxm::EwaldExclusionType::Table;
126 kernelSetup.ewaldExclusionType = options.useTabulatedEwaldCorr
127 ? Nbnxm::EwaldExclusionType::Table
128 : Nbnxm::EwaldExclusionType::Analytical;
134 void NbvSetupUtil::setParticleInfoAllVdv(const size_t numParticles)
137 particleInfoAllVdw_.resize(numParticles);
138 for (size_t particleI = 0; particleI < numParticles; particleI++)
140 particleInfoAllVdw_[particleI] |= gmx::sc_atomInfo_HasVdw;
141 particleInfoAllVdw_[particleI] |= gmx::sc_atomInfo_HasCharge;
145 void NbvSetupUtil::setNonBondedParameters(const std::vector<ParticleType>& particleTypes,
146 const NonBondedInteractionMap& nonBondedInteractionMap)
148 /* Todo: Refactor nbnxm to take nonbondedParameters_ directly
150 * initial self-handling of combination rules
151 * size: 2*(numParticleTypes^2)
153 nonbondedParameters_.reserve(2 * particleTypes.size() * particleTypes.size());
155 constexpr real c6factor = 6.0;
156 constexpr real c12factor = 12.0;
158 for (const ParticleType& particleType1 : particleTypes)
160 for (const ParticleType& particleType2 : particleTypes)
162 nonbondedParameters_.push_back(
163 nonBondedInteractionMap.getC6(particleType1.name(), particleType2.name()) * c6factor);
164 nonbondedParameters_.push_back(
165 nonBondedInteractionMap.getC12(particleType1.name(), particleType2.name()) * c12factor);
170 void NbvSetupUtil::setAtomProperties(const std::vector<int>& particleTypeIdOfAllParticles,
171 const std::vector<real>& charges)
173 gmxForceCalculator_->nbv_->setAtomProperties(particleTypeIdOfAllParticles, charges, particleInfoAllVdw_);
176 //! Sets up and returns a Nbnxm object for the given options and system
177 void NbvSetupUtil::setupNbnxmInstance(const size_t numParticleTypes, const NBKernelOptions& options)
179 const auto pinPolicy = (options.useGpu ? gmx::PinningPolicy::PinnedIfSupported
180 : gmx::PinningPolicy::CannotBePinned);
181 const int numThreads = options.numOpenMPThreads;
182 // Note: the options and Nbnxm combination rule enums values should match
183 const int combinationRule = static_cast<int>(options.ljCombinationRule);
185 checkKernelSetup(options); // throws exception is setup is invalid
187 Nbnxm::KernelSetup kernelSetup = getKernelSetup(options);
189 PairlistParams pairlistParams(kernelSetup.kernelType, false, options.pairlistCutoff, false);
190 Nbnxm::GridSet gridSet(
191 PbcType::Xyz, false, nullptr, nullptr, pairlistParams.pairlistType, false, numThreads, pinPolicy);
192 auto pairlistSets = std::make_unique<PairlistSets>(pairlistParams, false, 0);
193 auto pairSearch = std::make_unique<PairSearch>(
194 PbcType::Xyz, false, nullptr, nullptr, pairlistParams.pairlistType, false, numThreads, pinPolicy);
196 auto atomData = std::make_unique<nbnxn_atomdata_t>(pinPolicy,
198 kernelSetup.kernelType,
201 nonbondedParameters_,
205 // Put everything together
206 auto nbv = std::make_unique<nonbonded_verlet_t>(
207 std::move(pairlistSets), std::move(pairSearch), std::move(atomData), kernelSetup, nullptr, nullptr);
209 gmxForceCalculator_->nbv_ = std::move(nbv);
212 //! Computes the Ewald splitting coefficient for Coulomb
213 static real ewaldCoeff(const real ewald_rtol, const real pairlistCutoff)
215 return calc_ewaldcoeff_q(pairlistCutoff, ewald_rtol);
218 void NbvSetupUtil::setupStepWorkload(const NBKernelOptions& options)
220 gmxForceCalculator_->stepWork_->computeForces = true;
221 gmxForceCalculator_->stepWork_->computeNonbondedForces = true;
223 if (options.computeVirialAndEnergy)
225 gmxForceCalculator_->stepWork_->computeVirial = true;
226 gmxForceCalculator_->stepWork_->computeEnergy = true;
230 void NbvSetupUtil::setupInteractionConst(const NBKernelOptions& options)
232 gmxForceCalculator_->interactionConst_->vdwtype = VanDerWaalsType::Cut;
233 gmxForceCalculator_->interactionConst_->vdw_modifier = InteractionModifiers::PotShift;
234 gmxForceCalculator_->interactionConst_->rvdw = options.pairlistCutoff;
236 switch (options.coulombType)
238 case CoulombType::Pme:
239 gmxForceCalculator_->interactionConst_->eeltype = CoulombInteractionType::Pme;
241 case CoulombType::Cutoff:
242 gmxForceCalculator_->interactionConst_->eeltype = CoulombInteractionType::Cut;
244 case CoulombType::ReactionField:
245 gmxForceCalculator_->interactionConst_->eeltype = CoulombInteractionType::RF;
247 case CoulombType::Count: throw InputException("Unsupported electrostatic interaction");
249 gmxForceCalculator_->interactionConst_->coulomb_modifier = InteractionModifiers::PotShift;
250 gmxForceCalculator_->interactionConst_->rcoulomb = options.pairlistCutoff;
251 // Note: values correspond to ic->coulomb_modifier = InteractionModifiers::PotShift
252 gmxForceCalculator_->interactionConst_->dispersion_shift.cpot =
253 -1.0 / gmx::power6(gmxForceCalculator_->interactionConst_->rvdw);
254 gmxForceCalculator_->interactionConst_->repulsion_shift.cpot =
255 -1.0 / gmx::power12(gmxForceCalculator_->interactionConst_->rvdw);
257 // These are the initialized values but we leave them here so that later
258 // these can become options.
259 gmxForceCalculator_->interactionConst_->epsilon_r = 1.0;
260 gmxForceCalculator_->interactionConst_->reactionFieldPermitivity = 1.0;
262 /* Set the Coulomb energy conversion factor */
263 if (gmxForceCalculator_->interactionConst_->epsilon_r != 0)
265 gmxForceCalculator_->interactionConst_->epsfac =
266 ONE_4PI_EPS0 / gmxForceCalculator_->interactionConst_->epsilon_r;
270 /* eps = 0 is infinite dieletric: no Coulomb interactions */
271 gmxForceCalculator_->interactionConst_->epsfac = 0;
275 gmxForceCalculator_->interactionConst_->epsilon_r,
276 gmxForceCalculator_->interactionConst_->reactionFieldPermitivity,
277 gmxForceCalculator_->interactionConst_->rcoulomb,
278 &gmxForceCalculator_->interactionConst_->reactionFieldCoefficient,
279 &gmxForceCalculator_->interactionConst_->reactionFieldShift);
281 if (EEL_PME_EWALD(gmxForceCalculator_->interactionConst_->eeltype))
283 // Ewald coefficients, we ignore the potential shift
284 gmxForceCalculator_->interactionConst_->ewaldcoeff_q = ewaldCoeff(1e-5, options.pairlistCutoff);
285 if (gmxForceCalculator_->interactionConst_->ewaldcoeff_q <= 0)
287 throw InputException("Ewald coefficient should be > 0");
289 gmxForceCalculator_->interactionConst_->coulombEwaldTables =
290 std::make_unique<EwaldCorrectionTables>();
291 init_interaction_const_tables(nullptr, gmxForceCalculator_->interactionConst_.get(), 0, 0);
295 void NbvSetupUtil::setupForceRec(const matrix& box)
297 assert((gmxForceCalculator_->forcerec_ && "Forcerec not initialized"));
298 gmxForceCalculator_->forcerec_->nbfp = nonbondedParameters_;
299 gmxForceCalculator_->forcerec_->shift_vec.resize(numShiftVectors);
300 calc_shifts(box, gmxForceCalculator_->forcerec_->shift_vec);
303 void NbvSetupUtil::setParticlesOnGrid(const std::vector<Vec3>& coordinates, const Box& box)
305 gmxForceCalculator_->setParticlesOnGrid(particleInfoAllVdw_, coordinates, box);
308 void NbvSetupUtil::constructPairList(ExclusionLists<int> exclusionLists)
310 gmx::ListOfLists<int> exclusions(std::move(exclusionLists.ListRanges),
311 std::move(exclusionLists.ListElements));
312 gmxForceCalculator_->nbv_->constructPairlist(
313 gmx::InteractionLocality::Local, exclusions, 0, gmxForceCalculator_->nrnb_.get());
317 std::unique_ptr<GmxForceCalculator> GmxSetupDirector::setupGmxForceCalculator(const SimulationState& system,
318 const NBKernelOptions& options)
320 NbvSetupUtil nbvSetupUtil;
321 nbvSetupUtil.setExecutionContext(options);
322 nbvSetupUtil.setNonBondedParameters(system.topology().getParticleTypes(),
323 system.topology().getNonBondedInteractionMap());
324 nbvSetupUtil.setParticleInfoAllVdv(system.topology().numParticles());
326 nbvSetupUtil.setupInteractionConst(options);
327 nbvSetupUtil.setupStepWorkload(options);
328 nbvSetupUtil.setupNbnxmInstance(system.topology().getParticleTypes().size(), options);
329 nbvSetupUtil.setParticlesOnGrid(system.coordinates(), system.box());
330 nbvSetupUtil.constructPairList(system.topology().exclusionLists());
331 nbvSetupUtil.setAtomProperties(system.topology().getParticleTypeIdOfAllParticles(),
332 system.topology().getCharges());
333 nbvSetupUtil.setupForceRec(system.box().legacyMatrix());
335 return nbvSetupUtil.getGmxForceCalculator();