Refactor virtual site interface
[alexxy/gromacs.git] / src / gromacs / modularsimulator / forceelement.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019,2020, 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 /*! \libinternal \file
36  * \brief Declares the force element for the modular simulator
37  *
38  * This element calculates the forces, with or without shells or
39  * flexible constraints.
40  *
41  * \author Pascal Merz <pascal.merz@me.com>
42  * \ingroup module_modularsimulator
43  */
44
45 #ifndef GMX_MODULARSIMULATOR_FORCEELEMENT_H
46 #define GMX_MODULARSIMULATOR_FORCEELEMENT_H
47
48 #include <array>
49
50 #include "gromacs/domdec/dlbtiming.h"
51 #include "gromacs/mdtypes/md_enums.h"
52 #include "gromacs/utility/real.h"
53
54 #include "modularsimulatorinterfaces.h"
55 #include "topologyholder.h"
56
57 struct gmx_enfrot;
58 struct gmx_shellfc_t;
59 struct gmx_wallcycle;
60 struct pull_t;
61 struct t_fcdata;
62 struct t_nrnb;
63
64 namespace gmx
65 {
66 class Awh;
67 class EnergyElement;
68 class FreeEnergyPerturbationElement;
69 class ImdSession;
70 class MDAtoms;
71 class MdrunScheduleWorkload;
72 class StatePropagatorData;
73 class VirtualSitesHandler;
74
75 /*! \libinternal
76  * \ingroup module_modularsimulator
77  * \brief Force element
78  *
79  * The force element manages the call to either
80  * do_force(...) or relax_shell_flexcon(...)
81  */
82 class ForceElement final :
83     public ISimulatorElement,
84     public ITopologyHolderClient,
85     public INeighborSearchSignallerClient,
86     public IEnergySignallerClient
87 {
88 public:
89     //! Constructor
90     ForceElement(StatePropagatorData*           statePropagatorData,
91                  EnergyElement*                 energyElement,
92                  FreeEnergyPerturbationElement* freeEnergyPerturbationElement,
93                  bool                           isVerbose,
94                  bool                           isDynamicBox,
95                  FILE*                          fplog,
96                  const t_commrec*               cr,
97                  const t_inputrec*              inputrec,
98                  const MDAtoms*                 mdAtoms,
99                  t_nrnb*                        nrnb,
100                  t_forcerec*                    fr,
101                  t_fcdata*                      fcd,
102                  gmx_wallcycle*                 wcycle,
103                  MdrunScheduleWorkload*         runScheduleWork,
104                  VirtualSitesHandler*           vsite,
105                  ImdSession*                    imdSession,
106                  pull_t*                        pull_work,
107                  Constraints*                   constr,
108                  const gmx_mtop_t*              globalTopology,
109                  gmx_enfrot*                    enforcedRotation);
110
111     /*! \brief Register force calculation for step / time
112      *
113      * @param step                 The step number
114      * @param time                 The time
115      * @param registerRunFunction  Function allowing to register a run function
116      */
117     void scheduleTask(Step step, Time time, const RegisterRunFunctionPtr& registerRunFunction) override;
118
119     //! Check that we got the local topology
120     void elementSetup() override;
121     //! Print some final output
122     void elementTeardown() override;
123
124 private:
125     //! ITopologyHolderClient implementation
126     void setTopology(const gmx_localtop_t* top) override;
127     //! INeighborSearchSignallerClient implementation
128     SignallerCallbackPtr registerNSCallback() override;
129     //! IEnergySignallerClient implementation
130     SignallerCallbackPtr registerEnergyCallback(EnergySignallerEvent event) override;
131     //! The actual do_force call
132     template<bool doShellFC>
133     void run(Step step, Time time, unsigned int flags);
134
135     //! The shell / FC helper struct
136     gmx_shellfc_t* shellfc_;
137     //! Whether shells or flexible constraints are present
138     const bool doShellFC_;
139
140     //! The next NS step
141     Step nextNSStep_;
142     //! The next energy calculation step
143     Step nextEnergyCalculationStep_;
144     //! The next energy calculation step
145     Step nextVirialCalculationStep_;
146     //! The next free energy calculation step
147     Step nextFreeEnergyCalculationStep_;
148
149     //! Pointer to the micro state
150     StatePropagatorData* statePropagatorData_;
151     //! Pointer to the energy element
152     EnergyElement* energyElement_;
153     //! Pointer to the free energy perturbation element
154     FreeEnergyPerturbationElement* freeEnergyPerturbationElement_;
155
156     //! The local topology - updated by Topology via Client system
157     const gmx_localtop_t* localTopology_;
158
159     //! Whether we're having a dynamic box
160     const bool isDynamicBox_;
161     //! Whether we're being verbose
162     const bool isVerbose_;
163     //! The number of shell relaxation steps we did
164     Step nShellRelaxationSteps_;
165
166     //! DD / DLB helper object
167     const DDBalanceRegionHandler ddBalanceRegionHandler_;
168
169     /* \brief The FEP lambda vector
170      *
171      * Used if FEP is off, since do_force
172      * requires lambda to be allocated anyway
173      */
174     std::array<real, efptNR> lambda_;
175
176     // Access to ISimulator data
177     //! Handles logging.
178     FILE* fplog_;
179     //! Handles communication.
180     const t_commrec* cr_;
181     //! Contains user input mdp options.
182     const t_inputrec* inputrec_;
183     //! Atom parameters for this domain.
184     const MDAtoms* mdAtoms_;
185     //! Manages flop accounting.
186     t_nrnb* nrnb_;
187     //! Manages wall cycle accounting.
188     gmx_wallcycle* wcycle_;
189     //! Parameters for force calculations.
190     t_forcerec* fr_;
191     //! Handles virtual sites.
192     VirtualSitesHandler* vsite_;
193     //! The Interactive Molecular Dynamics session.
194     ImdSession* imdSession_;
195     //! The pull work object.
196     pull_t* pull_work_;
197     //! Helper struct for force calculations.
198     t_fcdata* fcd_;
199     //! Schedule of work for each MD step for this task.
200     MdrunScheduleWorkload* runScheduleWork_;
201     //! Handles constraints.
202     Constraints* constr_;
203     //! Handles enforced rotation.
204     gmx_enfrot* enforcedRotation_;
205 };
206
207 } // namespace gmx
208
209 #endif // GMX_MODULARSIMULATOR_FORCEELEMENT_H