2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 #ifndef GMX_MDLIB_COUPLING_H
39 #define GMX_MDLIB_COUPLING_H
46 #include "gromacs/math/vectypes.h"
47 #include "gromacs/mdtypes/md_enums.h"
48 #include "gromacs/utility/enumerationhelpers.h"
49 #include "gromacs/utility/real.h"
52 struct gmx_enerdata_t;
71 /* Update the size of per-atom arrays (e.g. after DD re-partitioning,
72 which might increase the number of home atoms). */
74 void update_tcouple(int64_t step,
75 const t_inputrec* inputrec,
77 gmx_ekindata_t* ekind,
78 const t_extmass* MassQ,
80 gmx::ArrayRef<const unsigned short> cTC);
82 /* Update Parrinello-Rahman, to be called before the coordinate update */
83 void update_pcouple_before_coordinates(FILE* fplog,
85 const t_inputrec* inputrec,
87 matrix parrinellorahmanMu,
91 /* Update the box, to be called after the coordinate update.
92 * For Berendsen P-coupling, also calculates the scaling factor
93 * and scales the coordinates.
94 * When the deform option is used, scales coordinates and box here.
96 void update_pcouple_after_coordinates(FILE* fplog,
98 const t_inputrec* inputrec,
100 gmx::ArrayRef<const unsigned short> cFREEZE,
101 const matrix pressure,
102 const matrix forceVirial,
103 const matrix constraintVirial,
104 matrix pressureCouplingMu,
107 gmx::BoxDeformation* boxDeformation,
108 bool scaleCoordinates);
110 /* Return TRUE if OK, FALSE in case of Shake Error */
112 extern bool update_randomize_velocities(const t_inputrec* ir,
116 gmx::ArrayRef<const unsigned short> cTC,
117 gmx::ArrayRef<const real> invMass,
118 gmx::ArrayRef<gmx::RVec> v,
119 const gmx::Update* upd,
120 const gmx::Constraints* constr);
122 void berendsen_tcoupl(const t_inputrec* ir,
123 gmx_ekindata_t* ekind,
125 std::vector<double>& therm_integral); //NOLINT(google-runtime-references)
127 void andersen_tcoupl(const t_inputrec* ir,
131 gmx::ArrayRef<const unsigned short> cTC,
132 gmx::ArrayRef<const real> invMass,
133 gmx::ArrayRef<gmx::RVec> v,
135 const std::vector<bool>& randomize,
136 gmx::ArrayRef<const real> boltzfac);
138 void nosehoover_tcoupl(const t_grpopts* opts,
139 const gmx_ekindata_t* ekind,
141 gmx::ArrayRef<double> xi,
142 gmx::ArrayRef<double> vxi,
143 const t_extmass* MassQ);
145 void trotter_update(const t_inputrec* ir,
147 gmx_ekindata_t* ekind,
148 const gmx_enerdata_t* enerd,
152 gmx::ArrayRef<const unsigned short> cTC,
153 gmx::ArrayRef<const real> invMass,
154 const t_extmass* MassQ,
155 gmx::ArrayRef<std::vector<int>> trotter_seqlist,
156 TrotterSequence trotter_seqno);
158 gmx::EnumerationArray<TrotterSequence, std::vector<int>>
159 init_npt_vars(const t_inputrec* ir, t_state* state, t_extmass* Mass, bool bTrotter);
161 real NPT_energy(const t_inputrec* ir, const t_state* state, const t_extmass* MassQ);
162 /* computes all the pressure/tempertature control energy terms to get a conserved energy */
164 void vrescale_tcoupl(const t_inputrec* ir,
166 gmx_ekindata_t* ekind,
168 gmx::ArrayRef<double> therm_integral);
169 /* Compute temperature scaling. For V-rescale it is done in update. */
171 void rescale_velocities(const gmx_ekindata_t* ekind,
172 gmx::ArrayRef<const unsigned short> cTC,
175 gmx::ArrayRef<gmx::RVec> v);
176 /* Rescale the velocities with the scaling factor in ekind */
179 * \brief Compute the new annealing temperature for a temperature group
181 * \param inputrec The input record
182 * \param temperatureGroup The temperature group
183 * \param time The current time
184 * \return The new reference temperature for the group
186 real computeAnnealingTargetTemperature(const t_inputrec& inputrec, int temperatureGroup, real time);
188 //! Check whether we do simulated annealing.
189 bool doSimulatedAnnealing(const t_inputrec* ir);
191 //! Initialize simulated annealing.
192 bool initSimulatedAnnealing(t_inputrec* ir, gmx::Update* upd);
194 // TODO: This is the only function in update.h altering the inputrec
195 void update_annealing_target_temp(t_inputrec* ir, real t, gmx::Update* upd);
196 /* Set reference temp for simulated annealing at time t*/
198 real calc_temp(real ekin, real nrdf);
199 /* Calculate the temperature */
201 real calc_pres(PbcType pbcType, int nwall, const matrix box, const tensor ekin, const tensor vir, tensor pres);
202 /* Calculate the pressure tensor, returns the scalar pressure.
203 * The unit of pressure is bar.
206 void parrinellorahman_pcoupl(FILE* fplog,
208 const t_inputrec* ir,
218 /*! \brief Calculate the pressure coupling scaling matrix
220 * Used by Berendsen and C-Rescale pressure coupling, this function
221 * computes the current value of the scaling matrix. The template
222 * parameter determines the pressure coupling algorithm.
224 template<PressureCoupling pressureCouplingType>
225 void pressureCouplingCalculateScalingMatrix(FILE* fplog,
227 const t_inputrec* ir,
231 const matrix force_vir,
232 const matrix constraint_vir,
234 double* baros_integral);
236 /*! \brief Scale the box and coordinates
238 * Used by Berendsen and C-Rescale pressure coupling, this function scales
239 * the box, the positions, and the velocities (C-Rescale only) according to
240 * the scaling matrix mu. The template parameter determines the pressure
241 * coupling algorithm.
243 template<PressureCoupling pressureCouplingType>
244 void pressureCouplingScaleBoxAndCoordinates(const t_inputrec* ir,
250 gmx::ArrayRef<gmx::RVec> x,
251 gmx::ArrayRef<gmx::RVec> v,
252 gmx::ArrayRef<const unsigned short> cFREEZE,
254 bool scaleCoordinates);
256 void pleaseCiteCouplingAlgorithms(FILE* fplog, const t_inputrec& ir);
258 /*! \brief Generate a new kinetic energy for the v-rescale thermostat
260 * Generates a new value for the kinetic energy, according to
261 * Bussi et al JCP (2007), Eq. (A7)
263 * This is used by update_tcoupl(), and by the VRescaleThermostat of the modular
265 * \todo Move this to the VRescaleThermostat once the modular simulator becomes
266 * the default code path.
268 * \param[in] kk present value of the kinetic energy of the atoms to be thermalized (in
270 * \param[in] sigma target average value of the kinetic energy (ndeg k_b T/2) (in
271 * the same units as kk)
272 * \param[in] ndeg number of degrees of freedom of the atoms to be thermalized
273 * \param[in] taut relaxation time of the thermostat, in units of 'how often this
275 * \param[in] step the time step this routine is called on
276 * \param[in] seed the random number generator seed
277 * \return the new kinetic energy
279 real vrescale_resamplekin(real kk, real sigma, real ndeg, real taut, int64_t step, int64_t seed);
281 #endif // GMX_MDLIB_COUPLING_H