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 \file
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
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 FreeEnergyPerturbationElement;
68 class MdrunScheduleWorkload;
69 class StatePropagatorData;
72 * \ingroup module_modularsimulator
73 * \brief Shell & flex constraints element
75 * The ShellFCElement manages the call to relax_shell_flexcon(...)
77 class ShellFCElement final :
78 public ISimulatorElement,
79 public ITopologyHolderClient,
80 public INeighborSearchSignallerClient,
81 public IEnergySignallerClient
85 ShellFCElement(StatePropagatorData* statePropagatorData,
86 EnergyElement* energyElement,
87 FreeEnergyPerturbationElement* freeEnergyPerturbationElement,
92 const t_inputrec* inputrec,
93 const MDAtoms* mdAtoms,
97 gmx_wallcycle* wcycle,
98 MdrunScheduleWorkload* runScheduleWork,
100 ImdSession* imdSession,
103 const gmx_mtop_t* globalTopology,
104 gmx_enfrot* enforcedRotation);
106 /*! \brief Register shell / flex constraint calculation for step / time
108 * @param step The step number
109 * @param time The time
110 * @param registerRunFunction Function allowing to register a run function
112 void scheduleTask(Step step, Time time, const RegisterRunFunctionPtr& registerRunFunction) override;
114 //! Check that we got the local topology
115 void elementSetup() override;
116 //! Print some final output
117 void elementTeardown() override;
119 //! Whether either shells or flexible constraints are used
120 static bool doShellsOrFlexConstraints(const gmx_mtop_t* mtop, int nflexcon);
123 //! ITopologyHolderClient implementation
124 void setTopology(const gmx_localtop_t* top) override;
125 //! INeighborSearchSignallerClient implementation
126 SignallerCallbackPtr registerNSCallback() override;
127 //! IEnergySignallerClient implementation
128 SignallerCallbackPtr registerEnergyCallback(EnergySignallerEvent event) override;
129 //! The actual do_force call
130 void run(Step step, Time time, unsigned int flags);
132 //! The shell / FC helper struct
133 gmx_shellfc_t* shellfc_;
137 //! The next energy calculation step
138 Step nextEnergyCalculationStep_;
139 //! The next energy calculation step
140 Step nextVirialCalculationStep_;
141 //! The next free energy calculation step
142 Step nextFreeEnergyCalculationStep_;
144 //! Pointer to the micro state
145 StatePropagatorData* statePropagatorData_;
146 //! Pointer to the energy element
147 EnergyElement* energyElement_;
148 //! The free energy perturbation element
149 FreeEnergyPerturbationElement* freeEnergyPerturbationElement_;
151 //! The local topology - updated by Topology via Client system
152 const gmx_localtop_t* localTopology_;
154 //! Whether we're having a dynamic box
155 const bool isDynamicBox_;
156 //! Whether we're being verbose
157 const bool isVerbose_;
158 //! The number of shell relaxation steps we did
161 //! DD / DLB helper object
162 const DDBalanceRegionHandler ddBalanceRegionHandler_;
164 /* \brief The FEP lambda vector
166 * Used if FEP is off, since do_force
167 * requires lambda to be allocated anyway
169 std::array<real, efptNR> lambda_;
171 // Access to ISimulator data
174 //! Handles communication.
175 const t_commrec* cr_;
176 //! Contains user input mdp options.
177 const t_inputrec* inputrec_;
178 //! Atom parameters for this domain.
179 const MDAtoms* mdAtoms_;
180 //! Manages flop accounting.
182 //! Manages wall cycle accounting.
183 gmx_wallcycle* wcycle_;
184 //! Parameters for force calculations.
186 //! Handles virtual sites.
188 //! The Interactive Molecular Dynamics session.
189 ImdSession* imdSession_;
190 //! The pull work object.
192 //! Helper struct for force calculations.
194 //! Schedule of work for each MD step for this task.
195 MdrunScheduleWorkload* runScheduleWork_;
196 //! Handles constraints.
197 Constraints* constr_;
198 //! Handles enforced rotation.
199 gmx_enfrot* enforcedRotation_;
204 #endif // GMX_MODULARSIMULATOR_SHELLFCELEMENT_H