41ca4234921b1806023a09567779e0e775da0f37
[alexxy/gromacs.git] / src / gromacs / modularsimulator / domdechelper.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019,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 Declares the domain decomposition helper for the modular simulator
37  *
38  * \author Pascal Merz <pascal.merz@me.com>
39  * \ingroup module_modularsimulator
40  *
41  * This header is only used within the modular simulator module
42  */
43
44 #ifndef GMX_MODULARSIMULATOR_DOMDECHELPER_H
45 #define GMX_MODULARSIMULATOR_DOMDECHELPER_H
46
47 #include "modularsimulatorinterfaces.h"
48
49 struct gmx_localtop_t;
50 struct gmx_wallcycle;
51 struct pull_t;
52 struct t_commrec;
53 struct t_forcerec;
54 struct t_inputrec;
55 struct t_nrnb;
56
57 namespace gmx
58 {
59 class Constraints;
60 class FreeEnergyPerturbationData;
61 class ImdSession;
62 class MDAtoms;
63 class MDLogger;
64 class StatePropagatorData;
65 class TopologyHolder;
66 class VirtualSitesHandler;
67
68 //! \addtogroup module_modularsimulator
69 //! \{
70
71 //! The function type allowing to request a check of the number of bonded interactions
72 typedef std::function<void()> CheckBondedInteractionsCallback;
73
74 /*! \internal
75  * \brief Infrastructure element responsible for domain decomposition
76  *
77  * This encapsulates the function call to domain decomposition, which is
78  * important for performance but outside of the current scope of the modular
79  * simulator project. This relies on legacy data structures for the state
80  * and the topology.
81  *
82  * This element does not implement the ISimulatorElement interface, as
83  * the Simulator is calling it explicitly between task queue population
84  * steps. This allows elements to be aware of any changes before
85  * deciding what functionality they need to run.
86  */
87 class DomDecHelper final : public INeighborSearchSignallerClient
88 {
89 public:
90     //! Constructor
91     DomDecHelper(bool                            isVerbose,
92                  int                             verbosePrintInterval,
93                  StatePropagatorData*            statePropagatorData,
94                  FreeEnergyPerturbationData*     freeEnergyPerturbationData,
95                  TopologyHolder*                 topologyHolder,
96                  CheckBondedInteractionsCallback checkBondedInteractionsCallback,
97                  int                             nstglobalcomm,
98                  FILE*                           fplog,
99                  t_commrec*                      cr,
100                  const MDLogger&                 mdlog,
101                  Constraints*                    constr,
102                  const t_inputrec*               inputrec,
103                  MDAtoms*                        mdAtoms,
104                  t_nrnb*                         nrnb,
105                  gmx_wallcycle*                  wcycle,
106                  t_forcerec*                     fr,
107                  VirtualSitesHandler*            vsite,
108                  ImdSession*                     imdSession,
109                  pull_t*                         pull_work);
110
111     /*! \brief Run domain decomposition
112      *
113      * Does domain decomposition partitioning at neighbor searching steps
114      *
115      * \param step  The step number
116      * \param time  The time
117      */
118     void run(Step step, Time time);
119
120     /*! \brief The initial domain decomposition partitioning
121      *
122      */
123     void setup();
124
125 private:
126     //! INeighborSearchSignallerClient implementation
127     std::optional<SignallerCallback> registerNSCallback() override;
128
129     //! The next NS step
130     Step nextNSStep_;
131     //! Whether we're being verbose
132     const bool isVerbose_;
133     //! If we're being verbose, how often?
134     const int verbosePrintInterval_;
135     //! The global communication frequency
136     const int nstglobalcomm_;
137
138     // TODO: Clarify relationship to data objects and find a more robust alternative to raw pointers (#3583)
139     //! Pointer to the micro state
140     StatePropagatorData* statePropagatorData_;
141     //! Pointer to the free energy data
142     FreeEnergyPerturbationData* freeEnergyPerturbationData_;
143     //! Pointer to the topology
144     TopologyHolder* topologyHolder_;
145     //! Pointer to the ComputeGlobalsHelper object - to ask for # of bonded interaction checking
146     CheckBondedInteractionsCallback checkBondedInteractionsCallback_;
147
148     //! Helper function unifying the DD partitioning calls in setup() and run()
149     void partitionSystem(bool                     verbose,
150                          bool                     isMasterState,
151                          int                      nstglobalcomm,
152                          gmx_wallcycle*           wcycle,
153                          std::unique_ptr<t_state> localState,
154                          t_state*                 globalState);
155
156     // Access to ISimulator data
157     //! Handles logging.
158     FILE* fplog_;
159     //! Handles communication.
160     t_commrec* cr_;
161     //! Handles logging.
162     const MDLogger& mdlog_;
163     //! Handles constraints.
164     Constraints* constr_;
165     //! Contains user input mdp options.
166     const t_inputrec* inputrec_;
167     //! Atom parameters for this domain.
168     MDAtoms* mdAtoms_;
169     //! Manages flop accounting.
170     t_nrnb* nrnb_;
171     //! Manages wall cycle accounting.
172     gmx_wallcycle* wcycle_;
173     //! Parameters for force calculations.
174     t_forcerec* fr_;
175     //! Handles virtual sites.
176     VirtualSitesHandler* vsite_;
177     //! The Interactive Molecular Dynamics session.
178     ImdSession* imdSession_;
179     //! The pull work object.
180     pull_t* pull_work_;
181 };
182
183 //! \}
184 } // namespace gmx
185
186 #endif // GMX_MODULARSIMULATOR_DOMDECHELPER_H