Add class ListedForces
[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, 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
74 #include "gromacs/math/vectypes.h"
75 #include "gromacs/topology/ifunc.h"
76 #include "gromacs/utility/arrayref.h"
77 #include "gromacs/utility/basedefinitions.h"
78 #include "gromacs/utility/classhelpers.h"
79
80 struct bonded_threading_t;
81 struct gmx_enerdata_t;
82 struct gmx_ffparams_t;
83 struct gmx_grppairener_t;
84 struct gmx_localtop_t;
85 struct gmx_multisim_t;
86 class history_t;
87 class InteractionDefinitions;
88 struct t_commrec;
89 struct t_fcdata;
90 struct t_forcerec;
91 struct t_lambda;
92 struct t_mdatoms;
93 struct t_nrnb;
94 class t_state;
95
96 namespace gmx
97 {
98 class ForceOutputs;
99 class StepWorkload;
100 template<typename>
101 class ArrayRef;
102 } // namespace gmx
103
104 //! Type of CPU function to compute a bonded interaction.
105 using BondedFunction = real (*)(int              nbonds,
106                                 const t_iatom    iatoms[],
107                                 const t_iparams  iparams[],
108                                 const rvec       x[],
109                                 rvec4            f[],
110                                 rvec             fshift[],
111                                 const t_pbc*     pbc,
112                                 real             lambda,
113                                 real*            dvdlambda,
114                                 const t_mdatoms* md,
115                                 t_fcdata*        fcd,
116                                 int*             ddgatindex);
117
118 //! Getter for finding a callable CPU function to compute an \c ftype interaction.
119 BondedFunction bondedFunction(int ftype);
120
121 /*! \libinternal
122  * \brief Class for calculating listed interactions, uses OpenMP parallelization
123  */
124 class ListedForces
125 {
126 public:
127     /*! \brief Constructor
128      *
129      * \param[in] numEnergyGroups  The number of energy groups, used for storage of pair energies
130      * \param[in] numThreads       The number of threads used for computed listed interactions
131      * \param[in] fplog            Log file for printing env.var. override, can be nullptr
132      */
133     ListedForces(int numEnergyGroups, int numThreads, FILE* fplog);
134
135     //! Destructor which is actually default but in the source file to hide implementation classes
136     ~ListedForces();
137
138     /*! \brief Copy the listed interactions from \p idef and set up the thread parallelization
139      *
140      * \param[in] idef           The idef with all listed interactions to be computed on this rank
141      * \param[in] numAtomsForce  Force are, potentially, computed for atoms 0 to \p numAtomsForce
142      * \param[in] useGpu         Whether a GPU is used to compute (part of) the listed interactions
143      */
144     void setup(const InteractionDefinitions& idef, int numAtomsForce, bool useGpu);
145
146     /*! \brief Do all aspects of energy and force calculations for mdrun
147      * on the set of listed interactions
148      *
149      * xWholeMolecules only needs to contain whole molecules when orientation
150      * restraints need to be computed and can be empty otherwise.
151      */
152     void calculate(struct gmx_wallcycle*          wcycle,
153                    const matrix                   box,
154                    const t_lambda*                fepvals,
155                    const t_commrec*               cr,
156                    const gmx_multisim_t*          ms,
157                    const rvec                     x[],
158                    gmx::ArrayRef<const gmx::RVec> xWholeMolecules,
159                    history_t*                     hist,
160                    gmx::ForceOutputs*             forceOutputs,
161                    const t_forcerec*              fr,
162                    const struct t_pbc*            pbc,
163                    gmx_enerdata_t*                enerd,
164                    t_nrnb*                        nrnb,
165                    const real*                    lambda,
166                    const t_mdatoms*               md,
167                    int*                           global_atom_index,
168                    const gmx::StepWorkload&       stepWork);
169
170     //! Returns whether bonded interactions are assigned to the CPU
171     bool haveCpuBondeds() const;
172
173     /*! \brief Returns whether listed forces are computed on the CPU
174      *
175      * NOTE: the current implementation returns true if there are position restraints
176      * or any bonded interactions computed on the CPU.
177      */
178     bool haveCpuListedForces() const;
179
180     //! Returns true if there are position, distance or orientation restraints
181     bool haveRestraints() const;
182
183     //! Returns a reference to the force calculation data
184     const t_fcdata& fcdata() const { return *fcdata_; }
185
186     //! Returns a reference to the force calculation data
187     t_fcdata& fcdata() { return *fcdata_; }
188
189 private:
190     //! The interaction definitions
191     InteractionDefinitions const* idef_ = nullptr;
192     //! Thread parallelization setup, unique_ptr to avoid declaring bonded_threading_t
193     std::unique_ptr<bonded_threading_t> threading_;
194     //! Data for bonded tables and NMR restraining, needs to be refactored
195     std::unique_ptr<t_fcdata> fcdata_;
196     //! Force buffer for free-energy forces
197     std::vector<real> forceBufferLambda_;
198     //! Shift force buffer for free-energy forces
199     std::vector<gmx::RVec> shiftForceBufferLambda_;
200
201     GMX_DISALLOW_COPY_MOVE_AND_ASSIGN(ListedForces);
202 };
203
204 #endif