160612f8e8228fd2d170d962b7924878074c2b25
[alexxy/gromacs.git] / api / nblib / nbnxmsetuphelpers.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 Utilities to setup GROMACS data structures for non-bonded 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 "gromacs/ewald/ewald_utils.h"
44 #include "gromacs/gpu_utils/device_stream_manager.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/gpu_data_mgmt.h"
54 #include "gromacs/nbnxm/nbnxm_gpu.h"
55 #include "gromacs/nbnxm/nbnxm.h"
56 #include "gromacs/nbnxm/nbnxm_simd.h"
57 #include "gromacs/nbnxm/pairlistset.h"
58 #include "gromacs/nbnxm/pairlistsets.h"
59 #include "gromacs/nbnxm/pairsearch.h"
60 #include "gromacs/pbcutil/pbc.h"
61 #include "gromacs/utility/logger.h"
62 #include "gromacs/utility/smalloc.h"
63 #include "nblib/exception.h"
64 #include "nblib/kerneloptions.h"
65 #include "nblib/particletype.h"
66
67 #include "nbnxmsetuphelpers.h"
68
69 namespace nblib
70 {
71
72 int64_t findNumEnergyGroups(gmx::ArrayRef<int64_t> particleInteractionFlags)
73 {
74     auto groupId = [](int code1, int code2) {
75         return (code1 & gmx::sc_atomInfo_EnergyGroupIdMask) < (code2 & gmx::sc_atomInfo_EnergyGroupIdMask);
76     };
77
78     int maxElement = *std::max_element(
79             std::begin(particleInteractionFlags), std::end(particleInteractionFlags), groupId);
80     return ((maxElement + 1) & gmx::sc_atomInfo_EnergyGroupIdMask);
81 }
82
83 Nbnxm::KernelType translateBenchmarkEnum(const SimdKernels& kernel)
84 {
85     int kernelInt = static_cast<int>(kernel);
86     return static_cast<Nbnxm::KernelType>(kernelInt);
87 }
88
89 void checkKernelSetupSimd(const SimdKernels nbnxmSimd)
90 {
91     if (nbnxmSimd >= SimdKernels::Count || nbnxmSimd == SimdKernels::SimdAuto)
92     {
93         throw InputException("Need a valid kernel SIMD type");
94     }
95     // Check SIMD support
96     if ((nbnxmSimd != SimdKernels::SimdNo && !GMX_SIMD)
97 #ifndef GMX_NBNXN_SIMD_4XN
98         || nbnxmSimd == SimdKernels::Simd4XM
99 #endif
100 #ifndef GMX_NBNXN_SIMD_2XNN
101         || nbnxmSimd == SimdKernels::Simd2XMM
102 #endif
103     )
104     {
105         throw InputException("The requested SIMD kernel was not set up at configuration time");
106     }
107 }
108
109 Nbnxm::KernelSetup createKernelSetupCPU(const SimdKernels nbnxmSimd, const bool useTabulatedEwaldCorr)
110 {
111     checkKernelSetupSimd(nbnxmSimd);
112
113     Nbnxm::KernelSetup kernelSetup;
114
115     // The int enum options.nbnxnSimd is set up to match Nbnxm::KernelType + 1
116     kernelSetup.kernelType = translateBenchmarkEnum(nbnxmSimd);
117
118     // The plain-C kernel does not support analytical ewald correction
119     if (kernelSetup.kernelType == Nbnxm::KernelType::Cpu4x4_PlainC)
120     {
121         kernelSetup.ewaldExclusionType = Nbnxm::EwaldExclusionType::Table;
122     }
123     else
124     {
125         kernelSetup.ewaldExclusionType = useTabulatedEwaldCorr ? Nbnxm::EwaldExclusionType::Table
126                                                                : Nbnxm::EwaldExclusionType::Analytical;
127     }
128
129     return kernelSetup;
130 }
131
132 std::vector<int64_t> createParticleInfoAllVdw(const size_t numParticles)
133 {
134     std::vector<int64_t> particleInfoAllVdw(numParticles);
135     for (size_t particleI = 0; particleI < numParticles; particleI++)
136     {
137         particleInfoAllVdw[particleI] |= gmx::sc_atomInfo_HasVdw;
138         particleInfoAllVdw[particleI] |= gmx::sc_atomInfo_HasCharge;
139     }
140     return particleInfoAllVdw;
141 }
142
143 std::vector<real> createNonBondedParameters(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     std::vector<real> nonbondedParameters;
152     nonbondedParameters.reserve(2 * particleTypes.size() * particleTypes.size());
153
154     constexpr real c6factor  = 6.0;
155     constexpr real c12factor = 12.0;
156
157     for (const ParticleType& particleType1 : particleTypes)
158     {
159         for (const ParticleType& particleType2 : particleTypes)
160         {
161             nonbondedParameters.push_back(
162                     nonBondedInteractionMap.getC6(particleType1.name(), particleType2.name()) * c6factor);
163             nonbondedParameters.push_back(
164                     nonBondedInteractionMap.getC12(particleType1.name(), particleType2.name()) * c12factor);
165         }
166     }
167     return nonbondedParameters;
168 }
169
170 gmx::StepWorkload createStepWorkload()
171 {
172     gmx::StepWorkload stepWorkload;
173     stepWorkload.computeForces          = true;
174     stepWorkload.computeNonbondedForces = true;
175     stepWorkload.useGpuFBufferOps       = false;
176     stepWorkload.useGpuXBufferOps       = false;
177
178     return stepWorkload;
179 }
180
181 real ewaldCoeff(const real ewald_rtol, const real pairlistCutoff)
182 {
183     return calc_ewaldcoeff_q(pairlistCutoff, ewald_rtol);
184 }
185
186 interaction_const_t createInteractionConst(const NBKernelOptions& options)
187 {
188     interaction_const_t interactionConst;
189     interactionConst.vdwtype      = VanDerWaalsType::Cut;
190     interactionConst.vdw_modifier = InteractionModifiers::PotShift;
191     interactionConst.rvdw         = options.pairlistCutoff;
192
193     switch (options.coulombType)
194     {
195         case CoulombType::Pme: interactionConst.eeltype = CoulombInteractionType::Pme; break;
196         case CoulombType::Cutoff: interactionConst.eeltype = CoulombInteractionType::Cut; break;
197         case CoulombType::ReactionField:
198             interactionConst.eeltype = CoulombInteractionType::RF;
199             break;
200         case CoulombType::Count: throw InputException("Unsupported electrostatic interaction");
201     }
202     interactionConst.coulomb_modifier = InteractionModifiers::PotShift;
203     interactionConst.rcoulomb         = options.pairlistCutoff;
204     // Note: values correspond to ic->coulomb_modifier = eintmodPOTSHIFT
205     interactionConst.dispersion_shift.cpot = -1.0 / gmx::power6(interactionConst.rvdw);
206     interactionConst.repulsion_shift.cpot  = -1.0 / gmx::power12(interactionConst.rvdw);
207
208     // These are the initialized values but we leave them here so that later
209     // these can become options.
210     interactionConst.epsilon_r                = 1.0;
211     interactionConst.reactionFieldPermitivity = 1.0;
212
213     /* Set the Coulomb energy conversion factor */
214     if (interactionConst.epsilon_r != 0)
215     {
216         interactionConst.epsfac = ONE_4PI_EPS0 / interactionConst.epsilon_r;
217     }
218     else
219     {
220         /* eps = 0 is infinite dieletric: no Coulomb interactions */
221         interactionConst.epsfac = 0;
222     }
223
224     calc_rffac(nullptr,
225                interactionConst.epsilon_r,
226                interactionConst.reactionFieldPermitivity,
227                interactionConst.rcoulomb,
228                &interactionConst.reactionFieldCoefficient,
229                &interactionConst.reactionFieldShift);
230
231
232     if (EEL_PME_EWALD(interactionConst.eeltype))
233     {
234         // Ewald coefficients, we ignore the potential shift
235         interactionConst.ewaldcoeff_q = ewaldCoeff(1e-5, options.pairlistCutoff);
236         if (interactionConst.ewaldcoeff_q <= 0)
237         {
238             throw InputException("Ewald coefficient should be > 0");
239         }
240         interactionConst.coulombEwaldTables = std::make_unique<EwaldCorrectionTables>();
241         init_interaction_const_tables(nullptr, &interactionConst, 0, 0);
242     }
243     return interactionConst;
244 }
245
246 std::unique_ptr<nonbonded_verlet_t> createNbnxmCPU(const size_t              numParticleTypes,
247                                                    const NBKernelOptions&    options,
248                                                    int                       numEnergyGroups,
249                                                    gmx::ArrayRef<const real> nonbondedParameters)
250 {
251     const auto pinPolicy  = gmx::PinningPolicy::CannotBePinned;
252     const int  numThreads = options.numOpenMPThreads;
253     // Note: the options and Nbnxm combination rule enums values should match
254     const int combinationRule = static_cast<int>(options.ljCombinationRule);
255
256     Nbnxm::KernelSetup kernelSetup =
257             createKernelSetupCPU(options.nbnxmSimd, options.useTabulatedEwaldCorr);
258
259     PairlistParams pairlistParams(kernelSetup.kernelType, false, options.pairlistCutoff, false);
260
261     auto pairlistSets = std::make_unique<PairlistSets>(pairlistParams, false, 0);
262     auto pairSearch   = std::make_unique<PairSearch>(
263             PbcType::Xyz, false, nullptr, nullptr, pairlistParams.pairlistType, false, numThreads, pinPolicy);
264
265     // Needs to be called with the number of unique ParticleTypes
266     auto atomData = std::make_unique<nbnxn_atomdata_t>(pinPolicy,
267                                                        gmx::MDLogger(),
268                                                        kernelSetup.kernelType,
269                                                        combinationRule,
270                                                        numParticleTypes,
271                                                        nonbondedParameters,
272                                                        numEnergyGroups,
273                                                        numThreads);
274
275     // Put everything together
276     auto nbv = std::make_unique<nonbonded_verlet_t>(
277             std::move(pairlistSets), std::move(pairSearch), std::move(atomData), kernelSetup, nullptr, nullptr);
278
279     return nbv;
280 }
281
282 void setGmxNonBondedNThreads(int numThreads)
283 {
284     gmx_omp_nthreads_set(ModuleMultiThread::Pairsearch, numThreads);
285     gmx_omp_nthreads_set(ModuleMultiThread::Nonbonded, numThreads);
286 }
287
288 void updateForcerec(t_forcerec* forcerec, const matrix& box)
289 {
290     assert(forcerec != nullptr && "Forcerec not initialized");
291     forcerec->shift_vec.resize(numShiftVectors);
292     calc_shifts(box, forcerec->shift_vec);
293 }
294
295
296 } // namespace nblib