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 v-rescale thermostat for the modular simulator
38 * \author Pascal Merz <pascal.merz@me.com>
39 * \ingroup module_modularsimulator
42 #ifndef GMX_MODULARSIMULATOR_VRESCALETHERMOSTAT_H
43 #define GMX_MODULARSIMULATOR_VRESCALETHERMOSTAT_H
45 #include "gromacs/utility/arrayref.h"
47 #include "energyelement.h"
48 #include "modularsimulatorinterfaces.h"
49 #include "propagator.h"
57 * \ingroup module_modularsimulator
58 * \brief Element implementing the v-rescale thermostat
60 * This element takes a callback to the propagator and updates the velocity
61 * scaling factor according to the v-rescale thermostat.
63 class VRescaleThermostat final : public ISimulatorElement, public ICheckpointHelperClient
67 VRescaleThermostat(int nstcouple,
71 int numTemperatureGroups,
72 double couplingTimeStep,
73 const real* referenceTemperature,
74 const real* couplingTime,
75 const real* numDegreesOfFreedom,
76 EnergyElement* energyElement,
77 ArrayRef<real> lambdaView,
78 PropagatorCallbackPtr propagatorCallback,
79 const t_state* globalState,
83 /*! \brief Register run function for step / time
85 * @param step The step number
86 * @param time The time
87 * @param registerRunFunction Function allowing to register a run function
89 void scheduleTask(Step step, Time time, const RegisterRunFunctionPtr& registerRunFunction) override;
91 //! No element setup needed
92 void elementSetup() override {}
93 //! No element teardown needed
94 void elementTeardown() override {}
96 //! Getter for the thermostatIntegral
97 const std::vector<double>& thermostatIntegral() const;
100 //! The frequency at which the thermostat is applied
101 const int nstcouple_;
102 //! If != 0, offset the step at which the thermostat is applied
104 //! Whether we're using full step kinetic energy
105 const bool useFullStepKE_;
109 //! The number of temperature groups
110 const int numTemperatureGroups_;
111 //! The coupling time step - simulation time step x nstcouple_
112 const double couplingTimeStep_;
113 //! Coupling temperature per group
114 const std::vector<real> referenceTemperature_;
115 //! Coupling time per group
116 const std::vector<real> couplingTime_;
117 //! Number of degrees of freedom per group
118 const std::vector<real> numDegreesOfFreedom_;
119 //! Work exerted by thermostat
120 std::vector<double> thermostatIntegral_;
122 //! Pointer to the energy element (for ekindata)
123 EnergyElement* energyElement_;
125 //! View on the scaling factor of the propagator
126 ArrayRef<real> lambda_;
127 //! Callback to let propagator know that we updated lambda
128 PropagatorCallbackPtr propagatorCallback_;
130 //! Set new lambda value (at T-coupling steps)
131 void setLambda(Step step);
133 //! ICheckpointHelperClient implementation
134 void writeCheckpoint(t_state* localState, t_state* globalState) override;
139 #endif // GMX_MODULARSIMULATOR_VRESCALETHERMOSTAT_H