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