6794bb47debd747241e8fc99336e816c9e238330
[alexxy/gromacs.git] / api / nblib / gmxbackenddata.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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 /*! \inpublicapi \file
36  * \brief
37  * Implements nblib simulation box
38  *
39  * \author Victor Holanda <victor.holanda@cscs.ch>
40  * \author Joe Jordan <ejjordan@kth.se>
41  * \author Prashanth Kanduri <kanduri@cscs.ch>
42  * \author Sebastian Keller <keller@cscs.ch>
43  */
44 #ifndef NBLIB_GMXBACKENDDATA_H
45 #define NBLIB_GMXBACKENDDATA_H
46
47 #include "gromacs/gmxlib/nrnb.h"
48 #include "gromacs/mdtypes/enerdata.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/nbnxm.h"
53 #include "gromacs/timing/wallcycle.h"
54 #include "gromacs/utility/listoflists.h"
55 #include "gromacs/utility/range.h"
56 #include "nblib/kerneloptions.h"
57 #include "nblib/nbnxmsetuphelpers.h"
58
59 namespace nblib
60 {
61
62 /*! \brief GROMACS non-bonded force calculation backend
63  *
64  * This class encapsulates the various GROMACS data structures and their interplay
65  * from the NBLIB user. The class is a private member of the ForceCalculator and
66  * is not intended for the public interface.
67  *
68  * Handles the task of storing the simulation problem description using the internal
69  * representation used within GROMACS. It currently supports short range non-bonded
70  * interactions (PP) on a single node CPU.
71  *
72  */
73 class GmxBackendData
74 {
75 public:
76     GmxBackendData() = default;
77     GmxBackendData(const NBKernelOptions& options,
78                    int                    numEnergyGroups,
79                    gmx::ArrayRef<int>     exclusionRanges,
80                    gmx::ArrayRef<int>     exclusionElements)
81     {
82         // Set hardware params from the execution context
83         setGmxNonBondedNThreads(options.numOpenMPThreads);
84
85         // Set interaction constants struct
86         interactionConst_ = createInteractionConst(options);
87
88         // Set up StepWorkload data
89         stepWork_ = createStepWorkload();
90
91         // Set up gmx_enerdata_t (holds energy information)
92         enerd_ = gmx_enerdata_t{ numEnergyGroups, 0 };
93
94         // Construct pair lists
95         std::vector<int> exclusionRanges_(exclusionRanges.begin(), exclusionRanges.end());
96         std::vector<int> exclusionElements_(exclusionElements.begin(), exclusionElements.end());
97         exclusions_ = gmx::ListOfLists<int>(std::move(exclusionRanges_), std::move(exclusionElements_));
98     }
99
100     //! exclusions in gmx format
101     gmx::ListOfLists<int> exclusions_;
102
103     //! Non-Bonded Verlet object for force calculation
104     std::unique_ptr<nonbonded_verlet_t> nbv_;
105
106     //! Only shift_vec is used
107     t_forcerec forcerec_;
108
109     //! Parameters for various interactions in the system
110     interaction_const_t interactionConst_;
111
112     //! Tasks to perform in an MD Step
113     gmx::StepWorkload stepWork_;
114
115     gmx::SimulationWorkload simulationWork_;
116
117     //! Energies of different interaction types; currently only needed as an argument for dispatchNonbondedKernel
118     gmx_enerdata_t enerd_{ 1, 0 };
119
120     //! Non-bonded flop counter; currently only needed as an argument for dispatchNonbondedKernel
121     t_nrnb nrnb_;
122
123     //! Keep track of whether updatePairlist has been called at least once
124     bool updatePairlistCalled{ false };
125 };
126
127 } // namespace nblib
128 #endif // NBLIB_GMXBACKENDDATA_H