8501ff062a6ce9c22bf7cf02679a88fd60bc443c
[alexxy/gromacs.git] / api / nblib / gmxsetup.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2020, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  * \brief Translation layer to GROMACS data structures for force calculations.
37  *
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>
42  */
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/smalloc.h"
61 #include "nblib/exception.h"
62 #include "nblib/kerneloptions.h"
63 #include "nblib/particletype.h"
64 #include "nblib/simulationstate.h"
65
66 namespace nblib
67 {
68
69 //! Helper to translate between the different enumeration values.
70 static Nbnxm::KernelType translateBenchmarkEnum(const SimdKernels& kernel)
71 {
72     int kernelInt = static_cast<int>(kernel);
73     return static_cast<Nbnxm::KernelType>(kernelInt);
74 }
75
76 /*! \brief Checks the kernel setup
77  *
78  * Returns an error string when the kernel is not available.
79  */
80 static void checkKernelSetup(const NBKernelOptions& options)
81 {
82     if (options.nbnxmSimd >= SimdKernels::Count || options.nbnxmSimd == SimdKernels::SimdAuto)
83     {
84         throw InputException("Need a valid kernel SIMD type");
85     }
86     // Check SIMD support
87     if ((options.nbnxmSimd != SimdKernels::SimdNo && !GMX_SIMD)
88 #ifndef GMX_NBNXN_SIMD_4XN
89         || options.nbnxmSimd == SimdKernels::Simd4XM
90 #endif
91 #ifndef GMX_NBNXN_SIMD_2XNN
92         || options.nbnxmSimd == SimdKernels::Simd2XMM
93 #endif
94     )
95     {
96         throw InputException("The requested SIMD kernel was not set up at configuration time");
97     }
98 }
99
100 NbvSetupUtil::NbvSetupUtil() : gmxForceCalculator_(std::make_unique<GmxForceCalculator>()) {}
101
102 void NbvSetupUtil::setExecutionContext(const NBKernelOptions& options)
103 {
104     // Todo: find a more general way to initialize hardware
105     gmx_omp_nthreads_set(emntPairsearch, options.numOpenMPThreads);
106     gmx_omp_nthreads_set(emntNonbonded, options.numOpenMPThreads);
107 }
108
109 Nbnxm::KernelSetup NbvSetupUtil::getKernelSetup(const NBKernelOptions& options)
110 {
111     checkKernelSetup(options);
112
113     Nbnxm::KernelSetup kernelSetup;
114
115     // The int enum options.nbnxnSimd is set up to match Nbnxm::KernelType + 1
116     kernelSetup.kernelType = translateBenchmarkEnum(options.nbnxmSimd);
117     // The plain-C kernel does not support analytical ewald correction
118     if (kernelSetup.kernelType == Nbnxm::KernelType::Cpu4x4_PlainC)
119     {
120         kernelSetup.ewaldExclusionType = Nbnxm::EwaldExclusionType::Table;
121     }
122     else
123     {
124         kernelSetup.ewaldExclusionType = options.useTabulatedEwaldCorr
125                                                  ? Nbnxm::EwaldExclusionType::Table
126                                                  : Nbnxm::EwaldExclusionType::Analytical;
127     }
128
129     return kernelSetup;
130 }
131
132 void NbvSetupUtil::setParticleInfoAllVdv(const size_t numParticles)
133
134 {
135     particleInfoAllVdw_.resize(numParticles);
136     for (size_t particleI = 0; particleI < numParticles; particleI++)
137     {
138         SET_CGINFO_HAS_VDW(particleInfoAllVdw_[particleI]);
139         SET_CGINFO_HAS_Q(particleInfoAllVdw_[particleI]);
140     }
141 }
142
143 void NbvSetupUtil::setNonBondedParameters(const std::vector<ParticleType>& particleTypes,
144                                           const NonBondedInteractionMap&   nonBondedInteractionMap)
145 {
146     /* Todo: Refactor nbnxm to take nonbondedParameters_ directly
147      *
148      * initial self-handling of combination rules
149      * size: 2*(numParticleTypes^2)
150      */
151     nonbondedParameters_.reserve(2 * particleTypes.size() * particleTypes.size());
152
153     constexpr real c6factor  = 6.0;
154     constexpr real c12factor = 12.0;
155
156     for (const ParticleType& particleType1 : particleTypes)
157     {
158         for (const ParticleType& particleType2 : particleTypes)
159         {
160             nonbondedParameters_.push_back(
161                     nonBondedInteractionMap.getC6(particleType1.name(), particleType2.name()) * c6factor);
162             nonbondedParameters_.push_back(
163                     nonBondedInteractionMap.getC12(particleType1.name(), particleType2.name()) * c12factor);
164         }
165     }
166 }
167
168 void NbvSetupUtil::setAtomProperties(const std::vector<int>&  particleTypeIdOfAllParticles,
169                                      const std::vector<real>& charges)
170 {
171     gmxForceCalculator_->nbv_->setAtomProperties(particleTypeIdOfAllParticles, charges, particleInfoAllVdw_);
172 }
173
174 //! Sets up and returns a Nbnxm object for the given options and system
175 void NbvSetupUtil::setupNbnxmInstance(const size_t numParticleTypes, const NBKernelOptions& options)
176 {
177     const auto pinPolicy  = (options.useGpu ? gmx::PinningPolicy::PinnedIfSupported
178                                            : gmx::PinningPolicy::CannotBePinned);
179     const int  numThreads = options.numOpenMPThreads;
180     // Note: the options and Nbnxm combination rule enums values should match
181     const int combinationRule = static_cast<int>(options.ljCombinationRule);
182
183     checkKernelSetup(options); // throws exception is setup is invalid
184
185     Nbnxm::KernelSetup kernelSetup = getKernelSetup(options);
186
187     PairlistParams pairlistParams(kernelSetup.kernelType, false, options.pairlistCutoff, false);
188     Nbnxm::GridSet gridSet(PbcType::Xyz, false, nullptr, nullptr, pairlistParams.pairlistType,
189                            false, numThreads, pinPolicy);
190     auto           pairlistSets = std::make_unique<PairlistSets>(pairlistParams, false, 0);
191     auto           pairSearch =
192             std::make_unique<PairSearch>(PbcType::Xyz, false, nullptr, nullptr,
193                                          pairlistParams.pairlistType, false, numThreads, pinPolicy);
194
195     auto atomData = std::make_unique<nbnxn_atomdata_t>(pinPolicy);
196
197     // Put everything together
198     auto nbv = std::make_unique<nonbonded_verlet_t>(std::move(pairlistSets), std::move(pairSearch),
199                                                     std::move(atomData), kernelSetup, nullptr,
200                                                     nullWallcycle);
201
202     // Needs to be called with the number of unique ParticleTypes
203     nbnxn_atomdata_init(gmx::MDLogger(), nbv->nbat.get(), kernelSetup.kernelType, combinationRule,
204                         numParticleTypes, nonbondedParameters_, 1, numThreads);
205
206     gmxForceCalculator_->nbv_ = std::move(nbv);
207 }
208
209 //! Computes the Ewald splitting coefficient for Coulomb
210 static real ewaldCoeff(const real ewald_rtol, const real pairlistCutoff)
211 {
212     return calc_ewaldcoeff_q(pairlistCutoff, ewald_rtol);
213 }
214
215 void NbvSetupUtil::setupStepWorkload(const NBKernelOptions& options)
216 {
217     gmxForceCalculator_->stepWork_->computeForces          = true;
218     gmxForceCalculator_->stepWork_->computeNonbondedForces = true;
219
220     if (options.computeVirialAndEnergy)
221     {
222         gmxForceCalculator_->stepWork_->computeVirial = true;
223         gmxForceCalculator_->stepWork_->computeEnergy = true;
224     }
225 }
226
227 void NbvSetupUtil::setupInteractionConst(const NBKernelOptions& options)
228 {
229     gmxForceCalculator_->interactionConst_->vdwtype      = evdwCUT;
230     gmxForceCalculator_->interactionConst_->vdw_modifier = eintmodPOTSHIFT;
231     gmxForceCalculator_->interactionConst_->rvdw         = options.pairlistCutoff;
232
233     switch (options.coulombType)
234     {
235         case CoulombType::Pme: gmxForceCalculator_->interactionConst_->eeltype = eelPME; break;
236         case CoulombType::Cutoff: gmxForceCalculator_->interactionConst_->eeltype = eelCUT; break;
237         case CoulombType::ReactionField:
238             gmxForceCalculator_->interactionConst_->eeltype = eelRF;
239             break;
240         case CoulombType::Count: throw InputException("Unsupported electrostatic interaction");
241     }
242     gmxForceCalculator_->interactionConst_->coulomb_modifier = eintmodPOTSHIFT;
243     gmxForceCalculator_->interactionConst_->rcoulomb         = options.pairlistCutoff;
244     // Note: values correspond to ic->coulomb_modifier = eintmodPOTSHIFT
245     gmxForceCalculator_->interactionConst_->dispersion_shift.cpot =
246             -1.0 / gmx::power6(gmxForceCalculator_->interactionConst_->rvdw);
247     gmxForceCalculator_->interactionConst_->repulsion_shift.cpot =
248             -1.0 / gmx::power12(gmxForceCalculator_->interactionConst_->rvdw);
249
250     // These are the initialized values but we leave them here so that later
251     // these can become options.
252     gmxForceCalculator_->interactionConst_->epsilon_r  = 1.0;
253     gmxForceCalculator_->interactionConst_->epsilon_rf = 1.0;
254
255     /* Set the Coulomb energy conversion factor */
256     if (gmxForceCalculator_->interactionConst_->epsilon_r != 0)
257     {
258         gmxForceCalculator_->interactionConst_->epsfac =
259                 ONE_4PI_EPS0 / gmxForceCalculator_->interactionConst_->epsilon_r;
260     }
261     else
262     {
263         /* eps = 0 is infinite dieletric: no Coulomb interactions */
264         gmxForceCalculator_->interactionConst_->epsfac = 0;
265     }
266
267     calc_rffac(nullptr, gmxForceCalculator_->interactionConst_->epsilon_r,
268                gmxForceCalculator_->interactionConst_->epsilon_rf,
269                gmxForceCalculator_->interactionConst_->rcoulomb,
270                &gmxForceCalculator_->interactionConst_->k_rf,
271                &gmxForceCalculator_->interactionConst_->c_rf);
272
273     if (EEL_PME_EWALD(gmxForceCalculator_->interactionConst_->eeltype))
274     {
275         // Ewald coefficients, we ignore the potential shift
276         gmxForceCalculator_->interactionConst_->ewaldcoeff_q = ewaldCoeff(1e-5, options.pairlistCutoff);
277         if (gmxForceCalculator_->interactionConst_->ewaldcoeff_q <= 0)
278         {
279             throw InputException("Ewald coefficient should be > 0");
280         }
281         gmxForceCalculator_->interactionConst_->coulombEwaldTables =
282                 std::make_unique<EwaldCorrectionTables>();
283         init_interaction_const_tables(nullptr, gmxForceCalculator_->interactionConst_.get(), 0);
284     }
285 }
286
287 void NbvSetupUtil::setupForceRec(const matrix& box)
288 {
289     assert((gmxForceCalculator_->forcerec_ && "Forcerec not initialized"));
290     gmxForceCalculator_->forcerec_->nbfp = nonbondedParameters_;
291     snew(gmxForceCalculator_->forcerec_->shift_vec, numShiftVectors);
292     calc_shifts(box, gmxForceCalculator_->forcerec_->shift_vec);
293 }
294
295 void NbvSetupUtil::setParticlesOnGrid(const std::vector<Vec3>& coordinates, const Box& box)
296 {
297     gmxForceCalculator_->setParticlesOnGrid(particleInfoAllVdw_, coordinates, box);
298 }
299
300 void NbvSetupUtil::constructPairList(const gmx::ListOfLists<int>& exclusions)
301 {
302     gmxForceCalculator_->nbv_->constructPairlist(gmx::InteractionLocality::Local, exclusions, 0,
303                                                  gmxForceCalculator_->nrnb_.get());
304 }
305
306
307 std::unique_ptr<GmxForceCalculator> GmxSetupDirector::setupGmxForceCalculator(const SimulationState& system,
308                                                                               const NBKernelOptions& options)
309 {
310     NbvSetupUtil nbvSetupUtil;
311     nbvSetupUtil.setExecutionContext(options);
312     nbvSetupUtil.setNonBondedParameters(system.topology().getParticleTypes(),
313                                         system.topology().getNonBondedInteractionMap());
314     nbvSetupUtil.setParticleInfoAllVdv(system.topology().numParticles());
315
316     nbvSetupUtil.setupInteractionConst(options);
317     nbvSetupUtil.setupStepWorkload(options);
318     nbvSetupUtil.setupNbnxmInstance(system.topology().getParticleTypes().size(), options);
319     nbvSetupUtil.setParticlesOnGrid(system.coordinates(), system.box());
320     nbvSetupUtil.constructPairList(system.topology().getGmxExclusions());
321     nbvSetupUtil.setAtomProperties(system.topology().getParticleTypeIdOfAllParticles(),
322                                    system.topology().getCharges());
323     nbvSetupUtil.setupForceRec(system.box().legacyMatrix());
324
325     return nbvSetupUtil.getGmxForceCalculator();
326 }
327
328 } // namespace nblib