Make temperature and pressure coupling enums enum classes
[alexxy/gromacs.git] / src / gromacs / modularsimulator / velocityscalingtemperaturecoupling.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 a velocity-scaling temperature coupling element for
37  * the modular simulator
38  *
39  * \author Pascal Merz <pascal.merz@me.com>
40  * \ingroup module_modularsimulator
41  *
42  * This header is only used within the modular simulator module
43  */
44
45 #ifndef GMX_MODULARSIMULATOR_VELOCITYSCALINGTEMPERATURECOUPLING_H
46 #define GMX_MODULARSIMULATOR_VELOCITYSCALINGTEMPERATURECOUPLING_H
47
48 #include "gromacs/utility/arrayref.h"
49
50 #include "energydata.h"
51 #include "modularsimulatorinterfaces.h"
52 #include "propagator.h"
53
54 struct t_commrec;
55
56 namespace gmx
57 {
58 class ITemperatureCouplingImpl;
59 class LegacySimulatorData;
60 struct TemperatureCouplingData;
61
62 //! Enum describing whether the thermostat is using full or half step kinetic energy
63 enum class UseFullStepKE
64 {
65     Yes,
66     No,
67     Count
68 };
69
70 //! Enum describing whether the thermostat is reporting conserved energy from the previous step
71 enum class ReportPreviousStepConservedEnergy
72 {
73     Yes,
74     No,
75     Count
76 };
77
78 /*! \internal
79  * \ingroup module_modularsimulator
80  * \brief Element implementing the a velocity-scaling thermostat
81  *
82  * This element takes a callback to the propagator and updates the velocity
83  * scaling factor according to the internal temperature coupling implementation.
84  *
85  * Note that the concrete implementation is handled by the concrete
86  * implementations of the ITemperatureCouplingImpl interface, while the element
87  * handles the scheduling and interfacing with other elements.
88  */
89 class VelocityScalingTemperatureCoupling final : public ISimulatorElement, public ICheckpointHelperClient
90 {
91 public:
92     //! Constructor
93     VelocityScalingTemperatureCoupling(int                               nstcouple,
94                                        int                               offset,
95                                        UseFullStepKE                     useFullStepKE,
96                                        ReportPreviousStepConservedEnergy reportPreviousConservedEnergy,
97                                        int64_t                           seed,
98                                        int                               numTemperatureGroups,
99                                        double                            couplingTimeStep,
100                                        const real*                       referenceTemperature,
101                                        const real*                       couplingTime,
102                                        const real*                       numDegreesOfFreedom,
103                                        EnergyData*                       energyData,
104                                        TemperatureCoupling               couplingType);
105
106     /*! \brief Register run function for step / time
107      *
108      * \param step                 The step number
109      * \param time                 The time
110      * \param registerRunFunction  Function allowing to register a run function
111      */
112     void scheduleTask(Step step, Time time, const RegisterRunFunction& registerRunFunction) override;
113
114     //! Sanity check at setup time
115     void elementSetup() override;
116     //! No element teardown needed
117     void elementTeardown() override {}
118
119     //! Contribution to the conserved energy (called by energy data)
120     [[nodiscard]] real conservedEnergyContribution() const;
121
122     //! Connect this to propagator
123     void connectWithPropagator(const PropagatorThermostatConnection& connectionData);
124
125     //! ICheckpointHelperClient write checkpoint implementation
126     void saveCheckpointState(std::optional<WriteCheckpointData> checkpointData, const t_commrec* cr) override;
127     //! ICheckpointHelperClient read checkpoint implementation
128     void restoreCheckpointState(std::optional<ReadCheckpointData> checkpointData, const t_commrec* cr) override;
129     //! ICheckpointHelperClient key implementation
130     const std::string& clientID() override;
131
132     /*! \brief Factory method implementation
133      *
134      * \param legacySimulatorData  Pointer allowing access to simulator level data
135      * \param builderHelper  ModularSimulatorAlgorithmBuilder helper object
136      * \param statePropagatorData  Pointer to the \c StatePropagatorData object
137      * \param energyData  Pointer to the \c EnergyData object
138      * \param freeEnergyPerturbationData  Pointer to the \c FreeEnergyPerturbationData object
139      * \param globalCommunicationHelper  Pointer to the \c GlobalCommunicationHelper object
140      * \param offset  The step offset at which the thermostat is applied
141      * \param useFullStepKE  Whether full step or half step KE is used
142      * \param reportPreviousStepConservedEnergy  Report the previous or the current step conserved energy
143      *
144      * \return  Pointer to the element to be added. Element needs to have been stored using \c storeElement
145      */
146     static ISimulatorElement*
147     getElementPointerImpl(LegacySimulatorData*                    legacySimulatorData,
148                           ModularSimulatorAlgorithmBuilderHelper* builderHelper,
149                           StatePropagatorData*                    statePropagatorData,
150                           EnergyData*                             energyData,
151                           FreeEnergyPerturbationData*             freeEnergyPerturbationData,
152                           GlobalCommunicationHelper*              globalCommunicationHelper,
153                           int                                     offset,
154                           UseFullStepKE                           useFullStepKE,
155                           ReportPreviousStepConservedEnergy reportPreviousStepConservedEnergy);
156
157 private:
158     //! The frequency at which the thermostat is applied
159     const int nstcouple_;
160     //! If != 0, offset the step at which the thermostat is applied
161     const int offset_;
162     //! Whether we're using full step kinetic energy
163     const UseFullStepKE useFullStepKE_;
164     //! Whether we are reporting the conserved energy from the previous step
165     const ReportPreviousStepConservedEnergy reportPreviousConservedEnergy_;
166
167     //! The number of temperature groups
168     const int numTemperatureGroups_;
169     //! The coupling time step - simulation time step x nstcouple_
170     const double couplingTimeStep_;
171     //! Coupling temperature per group
172     const std::vector<real> referenceTemperature_;
173     //! Coupling time per group
174     const std::vector<real> couplingTime_;
175     //! Number of degrees of freedom per group
176     const std::vector<real> numDegreesOfFreedom_;
177     //! Work exerted by thermostat per group
178     std::vector<double> temperatureCouplingIntegral_;
179     //! Work exerted by thermostat per group (backup from previous step)
180     std::vector<double> temperatureCouplingIntegralPreviousStep_;
181
182     // TODO: Clarify relationship to data objects and find a more robust alternative to raw pointers (#3583)
183     //! Pointer to the energy data (for ekindata)
184     EnergyData* energyData_;
185
186     //! Callback to let propagator know that we updated lambda
187     PropagatorCallback propagatorCallback_;
188
189     //! Set new lambda value (at T-coupling steps)
190     void setLambda(Step step);
191
192     //! The temperature coupling implementation
193     std::unique_ptr<ITemperatureCouplingImpl> temperatureCouplingImpl_;
194
195     //! CheckpointHelper identifier
196     const std::string identifier_ = "VelocityScalingTemperatureCoupling";
197     //! Helper function to read from / write to CheckpointData
198     template<CheckpointDataOperation operation>
199     void doCheckpointData(CheckpointData<operation>* checkpointData);
200 };
201
202 } // namespace gmx
203
204 #endif // GMX_MODULARSIMULATOR_VELOCITYSCALINGTEMPERATURECOUPLING_H