Add helper functions for setting up Nbnxm gpu object in nblib
[alexxy/gromacs.git] / api / nblib / nbnxmsetuphelpers.h
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 #ifndef NBLIB_NBNXMSETUPHELPERS_H
44 #define NBLIB_NBNXMSETUPHELPERS_H
45
46 #include "nblib/kerneloptions.h"
47 #include "nblib/interactions.h"
48 #include "nblib/particletype.h"
49
50 struct nonbonded_verlet_t;
51 struct t_forcerec;
52 struct t_nrnb;
53 struct interaction_const_t;
54 struct gmx_enerdata_t;
55 struct DeviceInformation;
56
57 namespace gmx
58 {
59 class StepWorkload;
60 class SimulationWorkload;
61 class DeviceStreamManager;
62 template<class T>
63 class ArrayRef;
64 } // namespace gmx
65
66 namespace Nbnxm
67 {
68 enum class KernelType;
69 struct KernelSetup;
70 } // namespace Nbnxm
71
72 namespace nblib
73 {
74
75 /*! \brief determine number of energy groups in an array of particle info flags
76  *
77  * Note: If the maximum energy group ID in the input is N, it is assumed that
78  * all the energy groups with IDs from 0...N-1 also exist.
79  */
80 int64_t findNumEnergyGroups(gmx::ArrayRef<int64_t> particleInteractionFlags);
81
82 //! Helper to translate between the different enumeration values.
83 Nbnxm::KernelType translateBenchmarkEnum(const SimdKernels& kernel);
84
85 /*! \brief Checks the kernel SIMD setup in CPU case
86  *
87  * Throws an exception when the kernel is not available.
88  */
89 void checkKernelSetupSimd(SimdKernels nbnxmSimd);
90
91 //! Creates and returns the kernel setup for CPU
92 Nbnxm::KernelSetup createKernelSetupCPU(const SimdKernels nbnxmSimd, const bool useTabulatedEwaldCorr);
93
94 //! Creates and returns the kernel setup for GPU
95 Nbnxm::KernelSetup createKernelSetupGPU(const bool useTabulatedEwaldCorr);
96
97 //! Create Particle info array to mark those that undergo VdV interaction
98 std::vector<int64_t> createParticleInfoAllVdw(size_t numParticles);
99
100 //! Create the non-bonded parameter vector in GROMACS format
101 std::vector<real> createNonBondedParameters(const std::vector<ParticleType>& particleTypes,
102                                             const NonBondedInteractionMap& nonBondedInteractionMap);
103
104 //! Create a step work object
105 gmx::StepWorkload createStepWorkload();
106
107 //! Create a SimulationWorkload object for use with createDeviceStreamManager
108 gmx::SimulationWorkload createSimulationWorkloadGpu();
109
110 //! Create a DeviceStreamManager; could be shared among multiple force calculators
111 std::shared_ptr<gmx::DeviceStreamManager> createDeviceStreamManager(const DeviceInformation& deviceInfo,
112                                                                     const gmx::SimulationWorkload& simulationWorkload);
113
114 //! Computes the Ewald splitting coefficient for Coulomb
115 real ewaldCoeff(real ewald_rtol, real pairlistCutoff);
116
117 //! Creates an interaction_const_t object from NBKernelOptions
118 interaction_const_t createInteractionConst(const NBKernelOptions& options);
119
120 //! Create nonbonded_verlet_t object
121 std::unique_ptr<nonbonded_verlet_t> createNbnxmCPU(size_t                    numParticleTypes,
122                                                    const NBKernelOptions&    options,
123                                                    int                       numEnergyGroups,
124                                                    gmx::ArrayRef<const real> nonbondedParameters);
125
126 //! Create nonbonded_verlet_gpu object
127 std::unique_ptr<nonbonded_verlet_t> createNbnxmGPU(size_t                     numParticleTypes,
128                                                    const NBKernelOptions&     options,
129                                                    const std::vector<real>&   nonbondedParameters,
130                                                    const interaction_const_t& interactionConst,
131                                                    const gmx::DeviceStreamManager& deviceStreamManager);
132
133 //! Set number of OpenMP threads in the GROMACS backend
134 void setGmxNonBondedNThreads(int numThreads);
135
136 //! Update the shift vectors in t_forcerec
137 void updateForcerec(t_forcerec* forcerec, const matrix& box);
138
139 } // namespace nblib
140
141 #endif // NBLIB_NBNXMSETUPHELPERS_H