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.
36 * \brief Declares the force element for the modular simulator
38 * \author Pascal Merz <pascal.merz@me.com>
39 * \ingroup module_modularsimulator
42 #ifndef GMX_MODULARSIMULATOR_FORCEELEMENT_H
43 #define GMX_MODULARSIMULATOR_FORCEELEMENT_H
47 #include "gromacs/domdec/dlbtiming.h"
48 #include "gromacs/mdtypes/md_enums.h"
49 #include "gromacs/utility/real.h"
51 #include "modularsimulatorinterfaces.h"
52 #include "topologyholder.h"
65 class MdrunScheduleWorkload;
66 class StatePropagatorData;
68 //! \addtogroup module_modularsimulator
72 * \brief Force element
74 * The force element manages the call to do_force(...)
76 class ForceElement final :
77 public ISimulatorElement,
78 public ITopologyHolderClient,
79 public INeighborSearchSignallerClient,
80 public IEnergySignallerClient
85 StatePropagatorData *statePropagatorData,
86 EnergyElement *energyElement,
90 const t_inputrec *inputrec,
91 const MDAtoms *mdAtoms,
95 gmx_wallcycle *wcycle,
96 MdrunScheduleWorkload *runScheduleWork,
98 ImdSession *imdSession,
101 /*! \brief Register force 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
108 Step step, Time time,
109 const RegisterRunFunctionPtr ®isterRunFunction) override;
111 //! Check that we got the local topology
112 void elementSetup() override;
113 //! No element teardown needed
114 void elementTeardown() override {}
117 //! ITopologyHolderClient implementation
118 void setTopology(const gmx_localtop_t *top) override;
119 //! INeighborSearchSignallerClient implementation
120 SignallerCallbackPtr registerNSCallback() override;
121 //! IEnergySignallerClient implementation
122 SignallerCallbackPtr registerEnergyCallback(EnergySignallerEvent event) override;
123 //! The actual do_force call
124 void run(Step step, Time time, unsigned int flags);
128 //! The next energy calculation step
129 Step nextEnergyCalculationStep_;
130 //! The next energy calculation step
131 Step nextVirialCalculationStep_;
132 //! The next free energy calculation step
133 Step nextFreeEnergyCalculationStep_;
135 //! Pointer to the micro state
136 StatePropagatorData *statePropagatorData_;
137 //! Pointer to the energy element
138 EnergyElement *energyElement_;
140 //! The local topology - updated by Topology via Client system
141 const gmx_localtop_t *localTopology_;
143 //! Whether we're having a dynamic box
144 const bool isDynamicBox_;
146 //! DD / DLB helper object
147 const DDBalanceRegionHandler ddBalanceRegionHandler_;
149 //! The FEP lambda vector (unused, but needs to be allocated)
150 std::array<real, efptNR> lambda_;
152 // Access to ISimulator data
155 //! Handles communication.
156 const t_commrec *cr_;
157 //! Contains user input mdp options.
158 const t_inputrec *inputrec_;
159 //! Atom parameters for this domain.
160 const MDAtoms *mdAtoms_;
161 //! Manages flop accounting.
163 //! Manages wall cycle accounting.
164 gmx_wallcycle *wcycle_;
165 //! Parameters for force calculations.
167 //! Handles virtual sites.
169 //! The Interactive Molecular Dynamics session.
170 ImdSession *imdSession_;
171 //! The pull work object.
173 //! Helper struct for force calculations.
175 //! Schedule of work for each MD step for this task.
176 MdrunScheduleWorkload *runScheduleWork_;
182 #endif // GMX_MODULARSIMULATOR_FORCEELEMENT_H