Remove preprocessor usage for atomInfo
[alexxy/gromacs.git] / api / nblib / gmxsetup.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
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/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"
67
68 namespace nblib
69 {
70
71 //! Helper to translate between the different enumeration values.
72 static Nbnxm::KernelType translateBenchmarkEnum(const SimdKernels& kernel)
73 {
74     int kernelInt = static_cast<int>(kernel);
75     return static_cast<Nbnxm::KernelType>(kernelInt);
76 }
77
78 /*! \brief Checks the kernel setup
79  *
80  * Returns an error string when the kernel is not available.
81  */
82 static void checkKernelSetup(const NBKernelOptions& options)
83 {
84     if (options.nbnxmSimd >= SimdKernels::Count || options.nbnxmSimd == SimdKernels::SimdAuto)
85     {
86         throw InputException("Need a valid kernel SIMD type");
87     }
88     // Check SIMD support
89     if ((options.nbnxmSimd != SimdKernels::SimdNo && !GMX_SIMD)
90 #ifndef GMX_NBNXN_SIMD_4XN
91         || options.nbnxmSimd == SimdKernels::Simd4XM
92 #endif
93 #ifndef GMX_NBNXN_SIMD_2XNN
94         || options.nbnxmSimd == SimdKernels::Simd2XMM
95 #endif
96     )
97     {
98         throw InputException("The requested SIMD kernel was not set up at configuration time");
99     }
100 }
101
102 NbvSetupUtil::NbvSetupUtil() : gmxForceCalculator_(std::make_unique<GmxForceCalculator>()) {}
103
104 void NbvSetupUtil::setExecutionContext(const NBKernelOptions& options)
105 {
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);
109 }
110
111 Nbnxm::KernelSetup NbvSetupUtil::getKernelSetup(const NBKernelOptions& options)
112 {
113     checkKernelSetup(options);
114
115     Nbnxm::KernelSetup kernelSetup;
116
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)
121     {
122         kernelSetup.ewaldExclusionType = Nbnxm::EwaldExclusionType::Table;
123     }
124     else
125     {
126         kernelSetup.ewaldExclusionType = options.useTabulatedEwaldCorr
127                                                  ? Nbnxm::EwaldExclusionType::Table
128                                                  : Nbnxm::EwaldExclusionType::Analytical;
129     }
130
131     return kernelSetup;
132 }
133
134 void NbvSetupUtil::setParticleInfoAllVdv(const size_t numParticles)
135
136 {
137     particleInfoAllVdw_.resize(numParticles);
138     for (size_t particleI = 0; particleI < numParticles; particleI++)
139     {
140         particleInfoAllVdw_[particleI] |= gmx::sc_atomInfo_HasVdw;
141         particleInfoAllVdw_[particleI] |= gmx::sc_atomInfo_HasCharge;
142     }
143 }
144
145 void NbvSetupUtil::setNonBondedParameters(const std::vector<ParticleType>& particleTypes,
146                                           const NonBondedInteractionMap&   nonBondedInteractionMap)
147 {
148     /* Todo: Refactor nbnxm to take nonbondedParameters_ directly
149      *
150      * initial self-handling of combination rules
151      * size: 2*(numParticleTypes^2)
152      */
153     nonbondedParameters_.reserve(2 * particleTypes.size() * particleTypes.size());
154
155     constexpr real c6factor  = 6.0;
156     constexpr real c12factor = 12.0;
157
158     for (const ParticleType& particleType1 : particleTypes)
159     {
160         for (const ParticleType& particleType2 : particleTypes)
161         {
162             nonbondedParameters_.push_back(
163                     nonBondedInteractionMap.getC6(particleType1.name(), particleType2.name()) * c6factor);
164             nonbondedParameters_.push_back(
165                     nonBondedInteractionMap.getC12(particleType1.name(), particleType2.name()) * c12factor);
166         }
167     }
168 }
169
170 void NbvSetupUtil::setAtomProperties(const std::vector<int>&  particleTypeIdOfAllParticles,
171                                      const std::vector<real>& charges)
172 {
173     gmxForceCalculator_->nbv_->setAtomProperties(particleTypeIdOfAllParticles, charges, particleInfoAllVdw_);
174 }
175
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)
178 {
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);
184
185     checkKernelSetup(options); // throws exception is setup is invalid
186
187     Nbnxm::KernelSetup kernelSetup = getKernelSetup(options);
188
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);
195
196     auto atomData = std::make_unique<nbnxn_atomdata_t>(pinPolicy,
197                                                        gmx::MDLogger(),
198                                                        kernelSetup.kernelType,
199                                                        combinationRule,
200                                                        numParticleTypes,
201                                                        nonbondedParameters_,
202                                                        1,
203                                                        numThreads);
204
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);
208
209     gmxForceCalculator_->nbv_ = std::move(nbv);
210 }
211
212 //! Computes the Ewald splitting coefficient for Coulomb
213 static real ewaldCoeff(const real ewald_rtol, const real pairlistCutoff)
214 {
215     return calc_ewaldcoeff_q(pairlistCutoff, ewald_rtol);
216 }
217
218 void NbvSetupUtil::setupStepWorkload(const NBKernelOptions& options)
219 {
220     gmxForceCalculator_->stepWork_->computeForces          = true;
221     gmxForceCalculator_->stepWork_->computeNonbondedForces = true;
222
223     if (options.computeVirialAndEnergy)
224     {
225         gmxForceCalculator_->stepWork_->computeVirial = true;
226         gmxForceCalculator_->stepWork_->computeEnergy = true;
227     }
228 }
229
230 void NbvSetupUtil::setupInteractionConst(const NBKernelOptions& options)
231 {
232     gmxForceCalculator_->interactionConst_->vdwtype      = VanDerWaalsType::Cut;
233     gmxForceCalculator_->interactionConst_->vdw_modifier = InteractionModifiers::PotShift;
234     gmxForceCalculator_->interactionConst_->rvdw         = options.pairlistCutoff;
235
236     switch (options.coulombType)
237     {
238         case CoulombType::Pme:
239             gmxForceCalculator_->interactionConst_->eeltype = CoulombInteractionType::Pme;
240             break;
241         case CoulombType::Cutoff:
242             gmxForceCalculator_->interactionConst_->eeltype = CoulombInteractionType::Cut;
243             break;
244         case CoulombType::ReactionField:
245             gmxForceCalculator_->interactionConst_->eeltype = CoulombInteractionType::RF;
246             break;
247         case CoulombType::Count: throw InputException("Unsupported electrostatic interaction");
248     }
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);
256
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;
261
262     /* Set the Coulomb energy conversion factor */
263     if (gmxForceCalculator_->interactionConst_->epsilon_r != 0)
264     {
265         gmxForceCalculator_->interactionConst_->epsfac =
266                 ONE_4PI_EPS0 / gmxForceCalculator_->interactionConst_->epsilon_r;
267     }
268     else
269     {
270         /* eps = 0 is infinite dieletric: no Coulomb interactions */
271         gmxForceCalculator_->interactionConst_->epsfac = 0;
272     }
273
274     calc_rffac(nullptr,
275                gmxForceCalculator_->interactionConst_->epsilon_r,
276                gmxForceCalculator_->interactionConst_->reactionFieldPermitivity,
277                gmxForceCalculator_->interactionConst_->rcoulomb,
278                &gmxForceCalculator_->interactionConst_->reactionFieldCoefficient,
279                &gmxForceCalculator_->interactionConst_->reactionFieldShift);
280
281     if (EEL_PME_EWALD(gmxForceCalculator_->interactionConst_->eeltype))
282     {
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)
286         {
287             throw InputException("Ewald coefficient should be > 0");
288         }
289         gmxForceCalculator_->interactionConst_->coulombEwaldTables =
290                 std::make_unique<EwaldCorrectionTables>();
291         init_interaction_const_tables(nullptr, gmxForceCalculator_->interactionConst_.get(), 0, 0);
292     }
293 }
294
295 void NbvSetupUtil::setupForceRec(const matrix& box)
296 {
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);
301 }
302
303 void NbvSetupUtil::setParticlesOnGrid(const std::vector<Vec3>& coordinates, const Box& box)
304 {
305     gmxForceCalculator_->setParticlesOnGrid(particleInfoAllVdw_, coordinates, box);
306 }
307
308 void NbvSetupUtil::constructPairList(ExclusionLists<int> exclusionLists)
309 {
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());
314 }
315
316
317 std::unique_ptr<GmxForceCalculator> GmxSetupDirector::setupGmxForceCalculator(const SimulationState& system,
318                                                                               const NBKernelOptions& options)
319 {
320     NbvSetupUtil nbvSetupUtil;
321     nbvSetupUtil.setExecutionContext(options);
322     nbvSetupUtil.setNonBondedParameters(system.topology().getParticleTypes(),
323                                         system.topology().getNonBondedInteractionMap());
324     nbvSetupUtil.setParticleInfoAllVdv(system.topology().numParticles());
325
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());
334
335     return nbvSetupUtil.getGmxForceCalculator();
336 }
337
338 } // namespace nblib