Refactor md_enums
[alexxy/gromacs.git] / src / gromacs / listed_forces / listed_forces.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2014,2015,2016,2017,2018 by the GROMACS development team.
5  * Copyright (c) 2019,2020,2021, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 /*! \defgroup module_listed_forces Interactions between lists of particles
37  * \ingroup group_mdrun
38  *
39  * \brief Handles computing energies and forces for listed
40  * interactions.
41  *
42  * Located here is the code for
43  * - computing energies and forces for interactions between a small
44      number of particles, e.g bonds, position restraints and listed
45      non-bonded interactions (e.g. 1-4).
46  * - high-level functions used by mdrun for computing a set of such
47      quantities
48  * - managing thread-wise decomposition, thread-local buffer output,
49      and reduction of output data across threads.
50  *
51  * \author Mark Abraham <mark.j.abraham@gmail.com>
52  * \author Berk Hess <hess@kth.se>
53  *
54  */
55 /*! \libinternal \file
56  *
57  * \brief This file contains declarations of high-level functions used
58  * by mdrun to compute energies and forces for listed interactions.
59  *
60  * Clients of libgromacs that want to evaluate listed interactions
61  * should call functions declared here.
62  *
63  * \author Mark Abraham <mark.j.abraham@gmail.com>
64  * \author Berk Hess <hess@kth.se>
65  *
66  * \inlibraryapi
67  * \ingroup module_listed_forces
68  */
69 #ifndef GMX_LISTED_FORCES_LISTED_FORCES_H
70 #define GMX_LISTED_FORCES_LISTED_FORCES_H
71
72 #include <memory>
73 #include <vector>
74
75 #include <bitset>
76
77 #include "gromacs/math/vectypes.h"
78 #include "gromacs/topology/idef.h"
79 #include "gromacs/topology/ifunc.h"
80 #include "gromacs/utility/basedefinitions.h"
81 #include "gromacs/utility/classhelpers.h"
82
83 struct bonded_threading_t;
84 struct gmx_enerdata_t;
85 struct gmx_ffparams_t;
86 struct gmx_grppairener_t;
87 struct gmx_localtop_t;
88 struct gmx_multisim_t;
89 class history_t;
90 struct t_commrec;
91 struct t_fcdata;
92 struct t_forcerec;
93 struct t_lambda;
94 struct t_mdatoms;
95 struct t_nrnb;
96 class t_state;
97 struct t_disresdata;
98 struct t_oriresdata;
99
100 namespace gmx
101 {
102 class ForceOutputs;
103 class StepWorkload;
104 template<typename>
105 class ArrayRef;
106 template<typename>
107 class ArrayRefWithPadding;
108 } // namespace gmx
109
110 //! Type of CPU function to compute a bonded interaction.
111 using BondedFunction = real (*)(int                 nbonds,
112                                 const t_iatom       iatoms[],
113                                 const t_iparams     iparams[],
114                                 const rvec          x[],
115                                 rvec4               f[],
116                                 rvec                fshift[],
117                                 const t_pbc*        pbc,
118                                 real                lambda,
119                                 gmx::ArrayRef<real> dvdlambda,
120                                 const t_mdatoms*    md,
121                                 t_fcdata*           fcd,
122                                 t_disresdata*       disresdata,
123                                 t_oriresdata*       oriresdata,
124                                 int*                ddgatindex);
125
126 //! Getter for finding a callable CPU function to compute an \c ftype interaction.
127 BondedFunction bondedFunction(int ftype);
128
129 /*! \libinternal
130  * \brief Class for calculating listed interactions, uses OpenMP parallelization
131  *
132  * Listed interactions can be divided over multiple instances of ListedForces
133  * using the selection flags passed to the constructor.
134  */
135 class ListedForces
136 {
137 public:
138     //! Enum for selecting groups of listed interaction types
139     enum class InteractionGroup : int
140     {
141         Pairs,     //!< Pair interactions
142         Dihedrals, //!< Dihedrals, including cmap
143         Angles,    //!< Angles
144         Rest,      //!< All listed interactions that are not any of the above
145         Count      //!< The number of items above
146     };
147
148     //! Type for specifying selections of groups of interaction types
149     using InteractionSelection = std::bitset<static_cast<int>(InteractionGroup::Count)>;
150
151     //! Returns a selection with all listed interaction types selected
152     static InteractionSelection interactionSelectionAll()
153     {
154         InteractionSelection is;
155         return is.flip();
156     }
157
158     /*! \brief Constructor
159      *
160      * \param[in] ffparams         The force field parameters
161      * \param[in] numEnergyGroups  The number of energy groups, used for storage of pair energies
162      * \param[in] numThreads       The number of threads used for computed listed interactions
163      * \param[in] interactionSelection  Select of interaction groups through bits set
164      * \param[in] fplog            Log file for printing env.var. override, can be nullptr
165      */
166     ListedForces(const gmx_ffparams_t& ffparams,
167                  int                   numEnergyGroups,
168                  int                   numThreads,
169                  InteractionSelection  interactionSelection,
170                  FILE*                 fplog);
171
172     //! Move constructor, default, but in the source file to hide implementation classes
173     ListedForces(ListedForces&& o) noexcept;
174
175     //! Destructor which is actually default but in the source file to hide implementation classes
176     ~ListedForces();
177
178     /*! \brief Copy the listed interactions from \p idef and set up the thread parallelization
179      *
180      * \param[in] domainIdef     Interaction definitions for all listed interactions to be computed on this domain/rank
181      * \param[in] numAtomsForce  Force are, potentially, computed for atoms 0 to \p numAtomsForce
182      * \param[in] useGpu         Whether a GPU is used to compute (part of) the listed interactions
183      */
184     void setup(const InteractionDefinitions& domainIdef, int numAtomsForce, bool useGpu);
185
186     /*! \brief Do all aspects of energy and force calculations for mdrun
187      * on the set of listed interactions
188      *
189      * xWholeMolecules only needs to contain whole molecules when orientation
190      * restraints need to be computed and can be empty otherwise.
191      */
192     void calculate(struct gmx_wallcycle*                     wcycle,
193                    const matrix                              box,
194                    const t_lambda*                           fepvals,
195                    const t_commrec*                          cr,
196                    const gmx_multisim_t*                     ms,
197                    gmx::ArrayRefWithPadding<const gmx::RVec> coordinates,
198                    gmx::ArrayRef<const gmx::RVec>            xWholeMolecules,
199                    t_fcdata*                                 fcdata,
200                    const history_t*                          hist,
201                    gmx::ForceOutputs*                        forceOutputs,
202                    const t_forcerec*                         fr,
203                    const struct t_pbc*                       pbc,
204                    gmx_enerdata_t*                           enerd,
205                    t_nrnb*                                   nrnb,
206                    gmx::ArrayRef<const real>                 lambda,
207                    const t_mdatoms*                          md,
208                    int*                                      global_atom_index,
209                    const gmx::StepWorkload&                  stepWork);
210
211     //! Returns whether bonded interactions are assigned to the CPU
212     bool haveCpuBondeds() const;
213
214     /*! \brief Returns whether listed forces are computed on the CPU
215      *
216      * NOTE: the current implementation returns true if there are position restraints
217      * or any bonded interactions computed on the CPU.
218      */
219     bool haveCpuListedForces(const t_fcdata& fcdata) const;
220
221     //! Returns true if there are position, distance or orientation restraints
222     bool haveRestraints(const t_fcdata& fcdata) const;
223
224 private:
225     //! Pointer to the interaction definitions
226     InteractionDefinitions const* idef_ = nullptr;
227     //! Interaction defintions used for storing selections
228     InteractionDefinitions idefSelection_;
229     //! Thread parallelization setup, unique_ptr to avoid declaring bonded_threading_t
230     std::unique_ptr<bonded_threading_t> threading_;
231     //! Tells which interactions to select for computation
232     const InteractionSelection interactionSelection_;
233     //! Force buffer for free-energy forces
234     std::vector<real> forceBufferLambda_;
235     //! Shift force buffer for free-energy forces
236     std::vector<gmx::RVec> shiftForceBufferLambda_;
237
238     GMX_DISALLOW_COPY_AND_ASSIGN(ListedForces);
239 };
240
241 #endif