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, 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_UPDATE_H
39 #define GMX_MDLIB_UPDATE_H
41 #include "gromacs/math/paddedvector.h"
42 #include "gromacs/math/vectypes.h"
43 #include "gromacs/mdtypes/md_enums.h"
44 #include "gromacs/timing/wallcycle.h"
45 #include "gromacs/utility/arrayref.h"
46 #include "gromacs/utility/basedefinitions.h"
47 #include "gromacs/utility/classhelpers.h"
48 #include "gromacs/utility/real.h"
51 struct gmx_ekindata_t;
52 struct gmx_enerdata_t;
63 /* Abstract type for update */
73 * \brief Contains data for update phase */
78 Update(const t_inputrec* ir, BoxDeformation* boxDeformation);
80 // TODO Get rid of getters when more free functions are incorporated as member methods
81 //! Returns handle to stochd_t struct
82 gmx_stochd_t* sd() const;
83 //! Returns pointer to PaddedVector xp
84 PaddedVector<gmx::RVec>* xp();
85 //! Returns handle to box deformation class
86 BoxDeformation* deform() const;
88 void setNumAtoms(int nAtoms);
91 //! Implementation type.
93 //! Implementation object.
94 PrivateImplPointer<Impl> impl_;
99 /* Update pre-computed constants that depend on the reference
100 * temperature for coupling.
102 * This could change e.g. in simulated annealing. */
103 void update_temperature_constants(gmx_stochd_t* sd, const t_inputrec* ir);
105 /* Update the size of per-atom arrays (e.g. after DD re-partitioning,
106 which might increase the number of home atoms). */
108 void update_tcouple(int64_t step,
109 const t_inputrec* inputrec,
111 gmx_ekindata_t* ekind,
112 const t_extmass* MassQ,
113 const t_mdatoms* md);
115 /* Update Parrinello-Rahman, to be called before the coordinate update */
116 void update_pcouple_before_coordinates(FILE* fplog,
118 const t_inputrec* inputrec,
120 matrix parrinellorahmanMu,
124 /* Update the box, to be called after the coordinate update.
125 * For Berendsen P-coupling, also calculates the scaling factor
126 * and scales the coordinates.
127 * When the deform option is used, scales coordinates and box here.
129 void update_pcouple_after_coordinates(FILE* fplog,
131 const t_inputrec* inputrec,
133 const matrix pressure,
134 const matrix forceVirial,
135 const matrix constraintVirial,
136 matrix pressureCouplingMu,
140 bool scaleCoordinates);
142 void update_coords(int64_t step,
143 const t_inputrec* inputrec, /* input record and box stuff */
146 gmx::ArrayRefWithPadding<const gmx::RVec> f, /* forces on home particles */
148 const gmx_ekindata_t* ekind,
152 const t_commrec* cr, /* these shouldn't be here -- need to think about it */
153 const gmx::Constraints* constr);
155 /* Return TRUE if OK, FALSE in case of Shake Error */
157 extern gmx_bool update_randomize_velocities(const t_inputrec* ir,
161 gmx::ArrayRef<gmx::RVec> v,
162 const gmx::Update* upd,
163 const gmx::Constraints* constr);
165 void update_sd_second_half(int64_t step,
166 real* dvdlambda, /* the contribution to be added to the bonded interactions */
167 const t_inputrec* inputrec, /* input record and box stuff */
172 gmx_wallcycle_t wcycle,
174 gmx::Constraints* constr,
178 void finish_update(const t_inputrec* inputrec,
181 gmx_wallcycle_t wcycle,
183 const gmx::Constraints* constr);
186 * Compute the partial kinetic energy for home particles;
187 * will be accumulated in the calling routine.
190 * Ekin = SUM(i) 0.5 m[i] v[i] (x) v[i]
192 * use v[i] = v[i] - u[i] when calculating temperature
194 * u must be accumulated already.
196 * Now also computes the contribution of the kinetic energy to the
202 void init_ekinstate(ekinstate_t* ekinstate, const t_inputrec* ir);
204 void update_ekinstate(ekinstate_t* ekinstate, const gmx_ekindata_t* ekind);
206 /*! \brief Restores data from \p ekinstate to \p ekind, then broadcasts it
207 to the rest of the simulation */
208 void restore_ekinstate_from_state(const t_commrec* cr, gmx_ekindata_t* ekind, const ekinstate_t* ekinstate);
210 void berendsen_tcoupl(const t_inputrec* ir,
211 gmx_ekindata_t* ekind,
213 std::vector<double>& therm_integral); //NOLINT(google-runtime-references)
215 void andersen_tcoupl(const t_inputrec* ir,
219 gmx::ArrayRef<gmx::RVec> v,
221 const std::vector<bool>& randomize,
222 gmx::ArrayRef<const real> boltzfac);
224 void nosehoover_tcoupl(const t_grpopts* opts,
225 const gmx_ekindata_t* ekind,
229 const t_extmass* MassQ);
231 void trotter_update(const t_inputrec* ir,
233 gmx_ekindata_t* ekind,
234 const gmx_enerdata_t* enerd,
238 const t_extmass* MassQ,
239 gmx::ArrayRef<std::vector<int>> trotter_seqlist,
242 std::array<std::vector<int>, ettTSEQMAX>
243 init_npt_vars(const t_inputrec* ir, t_state* state, t_extmass* Mass, gmx_bool bTrotter);
245 real NPT_energy(const t_inputrec* ir, const t_state* state, const t_extmass* MassQ);
246 /* computes all the pressure/tempertature control energy terms to get a conserved energy */
248 void vrescale_tcoupl(const t_inputrec* ir, int64_t step, gmx_ekindata_t* ekind, real dt, double therm_integral[]);
249 /* Compute temperature scaling. For V-rescale it is done in update. */
251 void rescale_velocities(const gmx_ekindata_t* ekind, const t_mdatoms* mdatoms, int start, int end, rvec v[]);
252 /* Rescale the velocities with the scaling factor in ekind */
254 //! Check whether we do simulated annealing.
255 bool doSimulatedAnnealing(const t_inputrec* ir);
257 //! Initialize simulated annealing.
258 bool initSimulatedAnnealing(t_inputrec* ir, gmx::Update* upd);
260 // TODO: This is the only function in update.h altering the inputrec
261 void update_annealing_target_temp(t_inputrec* ir, real t, gmx::Update* upd);
262 /* Set reference temp for simulated annealing at time t*/
264 real calc_temp(real ekin, real nrdf);
265 /* Calculate the temperature */
267 real calc_pres(PbcType pbcType, int nwall, const matrix box, const tensor ekin, const tensor vir, tensor pres);
268 /* Calculate the pressure tensor, returns the scalar pressure.
269 * The unit of pressure is bar.
272 void parrinellorahman_pcoupl(FILE* fplog,
274 const t_inputrec* ir,
282 gmx_bool bFirstStep);
284 void berendsen_pcoupl(FILE* fplog,
286 const t_inputrec* ir,
290 const matrix force_vir,
291 const matrix constraint_vir,
293 double* baros_integral);
295 void berendsen_pscale(const t_inputrec* ir,
302 const unsigned short cFREEZE[],
304 bool scaleCoordinates);
306 void pleaseCiteCouplingAlgorithms(FILE* fplog, const t_inputrec& ir);
308 /*! \brief Computes the atom range for a thread to operate on, ensuring SIMD aligned ranges
310 * \param[in] numThreads The number of threads to divide atoms over
311 * \param[in] threadIndex The thread to get the range for
312 * \param[in] numAtoms The total number of atoms (on this rank)
313 * \param[out] startAtom The start of the atom range
314 * \param[out] endAtom The end of the atom range, note that this is in general not a multiple of the SIMD width
316 void getThreadAtomRange(int numThreads, int threadIndex, int numAtoms, int* startAtom, int* endAtom);
318 /*! \brief Generate a new kinetic energy for the v-rescale thermostat
320 * Generates a new value for the kinetic energy, according to
321 * Bussi et al JCP (2007), Eq. (A7)
323 * This is used by update_tcoupl(), and by the VRescaleThermostat of the modular
325 * TODO: Move this to the VRescaleThermostat once the modular simulator becomes
326 * the default code path.
328 * @param kk present value of the kinetic energy of the atoms to be thermalized (in arbitrary units)
329 * @param sigma target average value of the kinetic energy (ndeg k_b T/2) (in the same units as kk)
330 * @param ndeg number of degrees of freedom of the atoms to be thermalized
331 * @param taut relaxation time of the thermostat, in units of 'how often this routine is called'
332 * @param step the time step this routine is called on
333 * @param seed the random number generator seed
334 * @return the new kinetic energy
336 real vrescale_resamplekin(real kk, real sigma, real ndeg, real taut, int64_t step, int64_t seed);