Domain decomposition and PME load balancing for modular simulator
[gromacs.git] / src / gromacs / modularsimulator / shellfcelement.h
blob81cae5617974cffe990f9277783b0afa02978f40
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019, 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.
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.
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.
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.
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.
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.
35 /*! \libinternal
36 * \brief Declares the shell / flex constraints element for the modular simulator
38 * \author Pascal Merz <pascal.merz@me.com>
39 * \ingroup module_modularsimulator
42 #ifndef GMX_MODULARSIMULATOR_SHELLFCELEMENT_H
43 #define GMX_MODULARSIMULATOR_SHELLFCELEMENT_H
45 #include "gromacs/domdec/dlbtiming.h"
47 #include "modularsimulatorinterfaces.h"
48 #include "topologyholder.h"
50 struct gmx_shellfc_t;
51 struct gmx_wallcycle;
52 struct pull_t;
53 struct t_fcdata;
54 struct t_nrnb;
56 namespace gmx
58 class Awh;
59 class EnergyElement;
60 class ImdSession;
61 class MDAtoms;
62 class MdScheduleWorkload;
63 class StatePropagatorData;
65 //! \addtogroup module_modularsimulator
66 //! \{
68 /*! \libinternal
69 * \brief Shell & flex constraints element
71 * The ShellFCElement manages the call to relax_shell_flexcon(...)
73 class ShellFCElement final :
74 public ISimulatorElement,
75 public ITopologyHolderClient,
76 public INeighborSearchSignallerClient,
77 public IEnergySignallerClient
79 public:
80 //! Constructor
81 ShellFCElement(
82 StatePropagatorData *statePropagatorData,
83 EnergyElement *energyElement,
84 bool isVerbose,
85 bool isDynamicBox,
86 FILE *fplog,
87 const t_commrec *cr,
88 const t_inputrec *inputrec,
89 const MDAtoms *mdAtoms,
90 t_nrnb *nrnb,
91 t_forcerec *fr,
92 t_fcdata *fcd,
93 gmx_wallcycle *wcycle,
94 MdScheduleWorkload *mdScheduleWork,
95 gmx_vsite_t *vsite,
96 ImdSession *imdSession,
97 pull_t *pull_work,
98 Constraints *constr,
99 const gmx_mtop_t *globalTopology);
101 /*! \brief Register shell / flex constraint calculation for step / time
103 * @param step The step number
104 * @param time The time
105 * @param registerRunFunction Function allowing to register a run function
107 void scheduleTask(
108 Step step, Time time,
109 const RegisterRunFunctionPtr &registerRunFunction) override;
111 //! Check that we got the local topology
112 void elementSetup() override;
113 //! Print some final output
114 void elementTeardown() override;
116 //! Whether either shells or flexible constraints are used
117 static bool doShellsOrFlexConstraints(const gmx_mtop_t *mtop, int nflexcon);
119 private:
120 //! ITopologyHolderClient implementation
121 void setTopology(const gmx_localtop_t *top) override;
122 //! INeighborSearchSignallerClient implementation
123 SignallerCallbackPtr registerNSCallback() override;
124 //! IEnergySignallerClient implementation
125 SignallerCallbackPtr registerEnergyCallback(EnergySignallerEvent event) override;
126 //! The actual do_force call
127 void run(Step step, Time time, unsigned int flags);
129 //! The shell / FC helper struct
130 gmx_shellfc_t *shellfc_;
132 //! The next NS step
133 Step nextNSStep_;
134 //! The next energy calculation step
135 Step nextEnergyCalculationStep_;
136 //! The next energy calculation step
137 Step nextVirialCalculationStep_;
138 //! The next free energy calculation step
139 Step nextFreeEnergyCalculationStep_;
141 //! Pointer to the micro state
142 StatePropagatorData *statePropagatorData_;
143 //! Pointer to the energy element
144 EnergyElement *energyElement_;
146 //! The local topology - updated by Topology via Client system
147 const gmx_localtop_t *localTopology_;
149 //! Whether we're having a dynamic box
150 const bool isDynamicBox_;
151 //! Whether we're being verbose
152 const bool isVerbose_;
153 //! The number of shell relaxation steps we did
154 Step nSteps_;
156 //! DD / DLB helper object
157 const DDBalanceRegionHandler ddBalanceRegionHandler_;
159 // Access to ISimulator data
160 //! Handles logging.
161 FILE *fplog_;
162 //! Handles communication.
163 const t_commrec *cr_;
164 //! Contains user input mdp options.
165 const t_inputrec *inputrec_;
166 //! Atom parameters for this domain.
167 const MDAtoms *mdAtoms_;
168 //! Manages flop accounting.
169 t_nrnb *nrnb_;
170 //! Manages wall cycle accounting.
171 gmx_wallcycle *wcycle_;
172 //! Parameters for force calculations.
173 t_forcerec *fr_;
174 //! Handles virtual sites.
175 gmx_vsite_t *vsite_;
176 //! The Interactive Molecular Dynamics session.
177 ImdSession *imdSession_;
178 //! The pull work object.
179 pull_t *pull_work_;
180 //! Helper struct for force calculations.
181 t_fcdata *fcd_;
182 //! Schedule of work for each MD step for this task.
183 MdScheduleWorkload *mdScheduleWork_;
184 //! Handles constraints.
185 Constraints *constr_;
188 //! \}
189 } // namespace gmx
191 #endif // GMX_MODULARSIMULATOR_SHELLFCELEMENT_H