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/real.h"
51 struct gmx_enerdata_t;
70 /* Update the size of per-atom arrays (e.g. after DD re-partitioning,
71 which might increase the number of home atoms). */
73 void update_tcouple(int64_t step,
74 const t_inputrec* inputrec,
76 gmx_ekindata_t* ekind,
77 const t_extmass* MassQ,
79 gmx::ArrayRef<const unsigned short> cTC);
81 /* Update Parrinello-Rahman, to be called before the coordinate update */
82 void update_pcouple_before_coordinates(FILE* fplog,
84 const t_inputrec* inputrec,
86 matrix parrinellorahmanMu,
90 /* Update the box, to be called after the coordinate update.
91 * For Berendsen P-coupling, also calculates the scaling factor
92 * and scales the coordinates.
93 * When the deform option is used, scales coordinates and box here.
95 void update_pcouple_after_coordinates(FILE* fplog,
97 const t_inputrec* inputrec,
99 gmx::ArrayRef<const unsigned short> cFREEZE,
100 const matrix pressure,
101 const matrix forceVirial,
102 const matrix constraintVirial,
103 matrix pressureCouplingMu,
106 gmx::BoxDeformation* boxDeformation,
107 bool scaleCoordinates);
109 /* Return TRUE if OK, FALSE in case of Shake Error */
111 extern bool update_randomize_velocities(const t_inputrec* ir,
115 gmx::ArrayRef<const unsigned short> cTC,
116 gmx::ArrayRef<const real> invMass,
117 gmx::ArrayRef<gmx::RVec> v,
118 const gmx::Update* upd,
119 const gmx::Constraints* constr);
121 void berendsen_tcoupl(const t_inputrec* ir,
122 gmx_ekindata_t* ekind,
124 std::vector<double>& therm_integral); //NOLINT(google-runtime-references)
126 void andersen_tcoupl(const t_inputrec* ir,
130 gmx::ArrayRef<const unsigned short> cTC,
131 gmx::ArrayRef<const real> invMass,
132 gmx::ArrayRef<gmx::RVec> v,
134 const std::vector<bool>& randomize,
135 gmx::ArrayRef<const real> boltzfac);
137 void nosehoover_tcoupl(const t_grpopts* opts,
138 const gmx_ekindata_t* ekind,
140 gmx::ArrayRef<double> xi,
141 gmx::ArrayRef<double> vxi,
142 const t_extmass* MassQ);
144 void trotter_update(const t_inputrec* ir,
146 gmx_ekindata_t* ekind,
147 const gmx_enerdata_t* enerd,
151 gmx::ArrayRef<const unsigned short> cTC,
152 gmx::ArrayRef<const real> invMass,
153 const t_extmass* MassQ,
154 gmx::ArrayRef<std::vector<int>> trotter_seqlist,
157 std::array<std::vector<int>, ettTSEQMAX>
158 init_npt_vars(const t_inputrec* ir, t_state* state, t_extmass* Mass, bool bTrotter);
160 real NPT_energy(const t_inputrec* ir, const t_state* state, const t_extmass* MassQ);
161 /* computes all the pressure/tempertature control energy terms to get a conserved energy */
163 void vrescale_tcoupl(const t_inputrec* ir,
165 gmx_ekindata_t* ekind,
167 gmx::ArrayRef<double> therm_integral);
168 /* Compute temperature scaling. For V-rescale it is done in update. */
170 void rescale_velocities(const gmx_ekindata_t* ekind,
171 gmx::ArrayRef<const unsigned short> cTC,
174 gmx::ArrayRef<gmx::RVec> v);
175 /* Rescale the velocities with the scaling factor in ekind */
177 //! Check whether we do simulated annealing.
178 bool doSimulatedAnnealing(const t_inputrec* ir);
180 //! Initialize simulated annealing.
181 bool initSimulatedAnnealing(t_inputrec* ir, gmx::Update* upd);
183 // TODO: This is the only function in update.h altering the inputrec
184 void update_annealing_target_temp(t_inputrec* ir, real t, gmx::Update* upd);
185 /* Set reference temp for simulated annealing at time t*/
187 real calc_temp(real ekin, real nrdf);
188 /* Calculate the temperature */
190 real calc_pres(PbcType pbcType, int nwall, const matrix box, const tensor ekin, const tensor vir, tensor pres);
191 /* Calculate the pressure tensor, returns the scalar pressure.
192 * The unit of pressure is bar.
195 void parrinellorahman_pcoupl(FILE* fplog,
197 const t_inputrec* ir,
207 /*! \brief Calculate the pressure coupling scaling matrix
209 * Used by Berendsen and C-Rescale pressure coupling, this function
210 * computes the current value of the scaling matrix. The template
211 * parameter determines the pressure coupling algorithm.
213 template<PressureCoupling pressureCouplingType>
214 void pressureCouplingCalculateScalingMatrix(FILE* fplog,
216 const t_inputrec* ir,
220 const matrix force_vir,
221 const matrix constraint_vir,
223 double* baros_integral);
225 /*! \brief Scale the box and coordinates
227 * Used by Berendsen and C-Rescale pressure coupling, this function scales
228 * the box, the positions, and the velocities (C-Rescale only) according to
229 * the scaling matrix mu. The template parameter determines the pressure
230 * coupling algorithm.
232 template<PressureCoupling pressureCouplingType>
233 void pressureCouplingScaleBoxAndCoordinates(const t_inputrec* ir,
239 gmx::ArrayRef<gmx::RVec> x,
240 gmx::ArrayRef<gmx::RVec> v,
241 gmx::ArrayRef<const unsigned short> cFREEZE,
243 bool scaleCoordinates);
245 void pleaseCiteCouplingAlgorithms(FILE* fplog, const t_inputrec& ir);
247 /*! \brief Generate a new kinetic energy for the v-rescale thermostat
249 * Generates a new value for the kinetic energy, according to
250 * Bussi et al JCP (2007), Eq. (A7)
252 * This is used by update_tcoupl(), and by the VRescaleThermostat of the modular
254 * \todo Move this to the VRescaleThermostat once the modular simulator becomes
255 * the default code path.
257 * \param[in] kk present value of the kinetic energy of the atoms to be thermalized (in
258 * arbitrary units) \param[in] sigma target average value of the kinetic energy (ndeg k_b T/2) (in
259 * the same units as kk) \param[in] ndeg number of degrees of freedom of the atoms to be
260 * thermalized \param[in] taut relaxation time of the thermostat, in units of 'how often this
261 * routine is called' \param[in] step the time step this routine is called on \param[in] seed the
262 * random number generator seed \return the new kinetic energy
264 real vrescale_resamplekin(real kk, real sigma, real ndeg, real taut, int64_t step, int64_t seed);
266 #endif // GMX_MDLIB_COUPLING_H