From: Joe Jordan Date: Fri, 11 Jun 2021 15:01:24 +0000 (+0000) Subject: Add createNbnxmCPU, usage, and test X-Git-Url: http://biod.pnpi.spb.ru/gitweb/?a=commitdiff_plain;h=c006875fb8407c3adc65babec07ee495ba60b42f;p=alexxy%2Fgromacs.git Add createNbnxmCPU, usage, and test --- diff --git a/api/nblib/gmxsetup.cpp b/api/nblib/gmxsetup.cpp index ba679f8451..27ead291e8 100644 --- a/api/nblib/gmxsetup.cpp +++ b/api/nblib/gmxsetup.cpp @@ -102,37 +102,8 @@ void NbvSetupUtil::setAtomProperties(const std::vector& particleTypeIdOfAl //! Sets up and returns a Nbnxm object for the given options and system void NbvSetupUtil::setupNbnxmInstance(const size_t numParticleTypes, const NBKernelOptions& options) { - const auto pinPolicy = (options.useGpu ? gmx::PinningPolicy::PinnedIfSupported - : gmx::PinningPolicy::CannotBePinned); - const int numThreads = options.numOpenMPThreads; - // Note: the options and Nbnxm combination rule enums values should match - const int combinationRule = static_cast(options.ljCombinationRule); - - checkKernelSetup(options.nbnxmSimd); // throws exception is setup is invalid - - Nbnxm::KernelSetup kernelSetup = getKernelSetup(options); - - PairlistParams pairlistParams(kernelSetup.kernelType, false, options.pairlistCutoff, false); - Nbnxm::GridSet gridSet( - PbcType::Xyz, false, nullptr, nullptr, pairlistParams.pairlistType, false, numThreads, pinPolicy); - auto pairlistSets = std::make_unique(pairlistParams, false, 0); - auto pairSearch = std::make_unique( - PbcType::Xyz, false, nullptr, nullptr, pairlistParams.pairlistType, false, numThreads, pinPolicy); - - auto atomData = std::make_unique(pinPolicy, - gmx::MDLogger(), - kernelSetup.kernelType, - combinationRule, - numParticleTypes, - nonbondedParameters_, - 1, - numThreads); - - // Put everything together - auto nbv = std::make_unique( - std::move(pairlistSets), std::move(pairSearch), std::move(atomData), kernelSetup, nullptr, nullptr); - - gmxForceCalculator_->nbv_ = std::move(nbv); + + gmxForceCalculator_->nbv_ = createNbnxmCPU(numParticleTypes, options, 1, nonbondedParameters_); } void NbvSetupUtil::setupStepWorkload(const NBKernelOptions& options) diff --git a/api/nblib/nbnxmsetuphelpers.cpp b/api/nblib/nbnxmsetuphelpers.cpp index 709ced3ce7..a1b4c9c22f 100644 --- a/api/nblib/nbnxmsetuphelpers.cpp +++ b/api/nblib/nbnxmsetuphelpers.cpp @@ -1,7 +1,7 @@ /* * This file is part of the GROMACS molecular simulation package. * - * Copyright (c) 2021, by the GROMACS development team, led by + * Copyright (c) 2020,2021, by the GROMACS development team, led by * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl, * and including many others, as listed in the AUTHORS file in the * top-level source directory and at http://www.gromacs.org. @@ -245,6 +245,41 @@ interaction_const_t createInteractionConst(const NBKernelOptions& options) return interactionConst; } +std::unique_ptr createNbnxmCPU(const size_t numParticleTypes, + const NBKernelOptions& options, + int numEnergyGroups, + gmx::ArrayRef nonbondedParameters) +{ + const auto pinPolicy = gmx::PinningPolicy::CannotBePinned; + const int numThreads = options.numOpenMPThreads; + // Note: the options and Nbnxm combination rule enums values should match + const int combinationRule = static_cast(options.ljCombinationRule); + + Nbnxm::KernelSetup kernelSetup = createKernelSetupCPU(options); + + PairlistParams pairlistParams(kernelSetup.kernelType, false, options.pairlistCutoff, false); + + auto pairlistSets = std::make_unique(pairlistParams, false, 0); + auto pairSearch = std::make_unique( + PbcType::Xyz, false, nullptr, nullptr, pairlistParams.pairlistType, false, numThreads, pinPolicy); + + // Needs to be called with the number of unique ParticleTypes + auto atomData = std::make_unique(pinPolicy, + gmx::MDLogger(), + kernelSetup.kernelType, + combinationRule, + numParticleTypes, + nonbondedParameters, + numEnergyGroups, + numThreads); + + // Put everything together + auto nbv = std::make_unique( + std::move(pairlistSets), std::move(pairSearch), std::move(atomData), kernelSetup, nullptr, nullptr); + + return nbv; +} + void setGmxNonBondedNThreads(int numThreads) { gmx_omp_nthreads_set(ModuleMultiThread::Pairsearch, numThreads); diff --git a/api/nblib/nbnxmsetuphelpers.h b/api/nblib/nbnxmsetuphelpers.h index decaa71e76..b3a2dc504b 100644 --- a/api/nblib/nbnxmsetuphelpers.h +++ b/api/nblib/nbnxmsetuphelpers.h @@ -107,6 +107,12 @@ real ewaldCoeff(real ewald_rtol, real pairlistCutoff); //! Creates an interaction_const_t object from NBKernelOptions interaction_const_t createInteractionConst(const NBKernelOptions& options); +//! Create nonbonded_verlet_t object +std::unique_ptr createNbnxmCPU(size_t numParticleTypes, + const NBKernelOptions& options, + int numEnergyGroups, + gmx::ArrayRef nonbondedParameters); + //! Set number of OpenMP threads in the GROMACS backend void setGmxNonBondedNThreads(int numThreads); diff --git a/api/nblib/tests/nbnxmsetup.cpp b/api/nblib/tests/nbnxmsetup.cpp index 5b02de6f00..a928ed4aa6 100644 --- a/api/nblib/tests/nbnxmsetup.cpp +++ b/api/nblib/tests/nbnxmsetup.cpp @@ -58,6 +58,7 @@ namespace test { namespace { + TEST(NbnxmSetupTest, findNumEnergyGroups) { std::vector v(10); @@ -178,6 +179,16 @@ TEST(NbnxmSetupTest, cannotCreateKernelSetupCPU4XM) } #endif +TEST(NbnxmSetupTest, CanCreateNbnxmCPU) +{ + size_t numParticles = 1; + NBKernelOptions nbKernelOptions; + nbKernelOptions.nbnxmSimd = SimdKernels::SimdNo; + int numEnergyGroups = 1; + std::vector nonbondedParameters = { 1, 1 }; + EXPECT_NO_THROW(createNbnxmCPU(numParticles, nbKernelOptions, numEnergyGroups, nonbondedParameters)); +} + } // namespace } // namespace test } // namespace nblib