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/nbnxmsetuphelpers.h"
66 #include "nblib/particletype.h"
67 #include "nblib/simulationstate.h"
72 NbvSetupUtil::NbvSetupUtil() : gmxForceCalculator_(std::make_unique<GmxForceCalculator>()) {}
74 void NbvSetupUtil::setExecutionContext(const NBKernelOptions& options)
76 setGmxNonBondedNThreads(options.numOpenMPThreads);
79 Nbnxm::KernelSetup NbvSetupUtil::getKernelSetup(const NBKernelOptions& options)
81 return createKernelSetupCPU(options);
84 void NbvSetupUtil::setParticleInfoAllVdv(const size_t numParticles)
87 particleInfoAllVdw_ = createParticleInfoAllVdv(numParticles);
90 void NbvSetupUtil::setNonBondedParameters(const std::vector<ParticleType>& particleTypes,
91 const NonBondedInteractionMap& nonBondedInteractionMap)
93 nonbondedParameters_ = createNonBondedParameters(particleTypes, nonBondedInteractionMap);
96 void NbvSetupUtil::setAtomProperties(const std::vector<int>& particleTypeIdOfAllParticles,
97 const std::vector<real>& charges)
99 gmxForceCalculator_->nbv_->setAtomProperties(particleTypeIdOfAllParticles, charges, particleInfoAllVdw_);
102 //! Sets up and returns a Nbnxm object for the given options and system
103 void NbvSetupUtil::setupNbnxmInstance(const size_t numParticleTypes, const NBKernelOptions& options)
106 gmxForceCalculator_->nbv_ = createNbnxmCPU(numParticleTypes, options, 1, nonbondedParameters_);
109 void NbvSetupUtil::setupStepWorkload(const NBKernelOptions& options)
111 gmxForceCalculator_->stepWork_ = std::make_unique<gmx::StepWorkload>(createStepWorkload(options));
114 void NbvSetupUtil::setupInteractionConst(const NBKernelOptions& options)
116 gmxForceCalculator_->interactionConst_ =
117 std::make_unique<interaction_const_t>(createInteractionConst(options));
120 void NbvSetupUtil::setupForceRec(const matrix& box)
122 updateForcerec(gmxForceCalculator_->forcerec_.get(), box);
125 void NbvSetupUtil::setParticlesOnGrid(const std::vector<Vec3>& coordinates, const Box& box)
127 gmxForceCalculator_->setParticlesOnGrid(particleInfoAllVdw_, coordinates, box);
130 void NbvSetupUtil::constructPairList(ExclusionLists<int> exclusionLists)
132 gmx::ListOfLists<int> exclusions(std::move(exclusionLists.ListRanges),
133 std::move(exclusionLists.ListElements));
134 gmxForceCalculator_->nbv_->constructPairlist(
135 gmx::InteractionLocality::Local, exclusions, 0, gmxForceCalculator_->nrnb_.get());
139 std::unique_ptr<GmxForceCalculator> GmxSetupDirector::setupGmxForceCalculator(const SimulationState& system,
140 const NBKernelOptions& options)
142 NbvSetupUtil nbvSetupUtil;
143 nbvSetupUtil.setExecutionContext(options);
144 nbvSetupUtil.setNonBondedParameters(system.topology().getParticleTypes(),
145 system.topology().getNonBondedInteractionMap());
146 nbvSetupUtil.setParticleInfoAllVdv(system.topology().numParticles());
148 nbvSetupUtil.setupInteractionConst(options);
149 nbvSetupUtil.setupStepWorkload(options);
150 nbvSetupUtil.setupNbnxmInstance(system.topology().getParticleTypes().size(), options);
151 nbvSetupUtil.setParticlesOnGrid(system.coordinates(), system.box());
152 nbvSetupUtil.constructPairList(system.topology().exclusionLists());
153 nbvSetupUtil.setAtomProperties(system.topology().getParticleTypeIdOfAllParticles(),
154 system.topology().getCharges());
155 nbvSetupUtil.setupForceRec(system.box().legacyMatrix());
157 return nbvSetupUtil.getGmxForceCalculator();