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,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 #ifndef GMX_MDLIB_UPDATE_H
38 #define GMX_MDLIB_UPDATE_H
40 #include "gromacs/math/paddedvector.h"
41 #include "gromacs/math/vectypes.h"
42 #include "gromacs/mdtypes/md_enums.h"
43 #include "gromacs/timing/wallcycle.h"
44 #include "gromacs/utility/arrayref.h"
45 #include "gromacs/utility/basedefinitions.h"
46 #include "gromacs/utility/classhelpers.h"
47 #include "gromacs/utility/real.h"
50 struct gmx_ekindata_t;
51 struct gmx_enerdata_t;
61 /* Abstract type for update */
71 * \brief Contains data for update phase */
76 Update(const t_inputrec *ir, BoxDeformation *boxDeformation);
78 // TODO Get rid of getters when more free functions are incorporated as member methods
79 //! Returns handle to stochd_t struct
80 gmx_stochd_t * sd() const;
81 //! Returns pointer to PaddedVector xp
82 PaddedVector<gmx::RVec> * xp();
83 //! Returns handle to box deformation class
84 BoxDeformation * deform() const;
86 void setNumAtoms(int nAtoms);
88 //! Implementation type.
90 //! Implementation object.
91 PrivateImplPointer<Impl> impl_;
96 /* Update pre-computed constants that depend on the reference
97 * temperature for coupling.
99 * This could change e.g. in simulated annealing. */
100 void update_temperature_constants(gmx_stochd_t *sd, const t_inputrec *ir);
102 /* Update the size of per-atom arrays (e.g. after DD re-partitioning,
103 which might increase the number of home atoms). */
105 void update_tcouple(int64_t step,
106 const t_inputrec *inputrec,
108 gmx_ekindata_t *ekind,
109 const t_extmass *MassQ,
113 /* Update Parrinello-Rahman, to be called before the coordinate update */
114 void update_pcouple_before_coordinates(FILE *fplog,
116 const t_inputrec *inputrec,
118 matrix parrinellorahmanMu,
122 /* Update the box, to be called after the coordinate update.
123 * For Berendsen P-coupling, also calculates the scaling factor
124 * and scales the coordinates.
125 * When the deform option is used, scales coordinates and box here.
127 void update_pcouple_after_coordinates(FILE *fplog,
129 const t_inputrec *inputrec,
131 const matrix pressure,
132 const matrix forceVirial,
133 const matrix constraintVirial,
134 const matrix parrinellorahmanMu,
139 void update_coords(int64_t step,
140 const t_inputrec *inputrec, /* input record and box stuff */
143 gmx::ArrayRefWithPadding<const gmx::RVec> f, /* forces on home particles */
145 const gmx_ekindata_t *ekind,
149 const t_commrec *cr, /* these shouldn't be here -- need to think about it */
150 const gmx::Constraints *constr);
152 /* Return TRUE if OK, FALSE in case of Shake Error */
154 extern gmx_bool update_randomize_velocities(const t_inputrec *ir, int64_t step, const t_commrec *cr,
156 gmx::ArrayRef<gmx::RVec> v,
157 const gmx::Update *upd,
158 const gmx::Constraints *constr);
160 void constrain_velocities(int64_t step,
161 real *dvdlambda, /* the contribution to be added to the bonded interactions */
164 gmx::Constraints *constr,
169 void constrain_coordinates(int64_t step,
170 real *dvdlambda, /* the contribution to be added to the bonded interactions */
174 gmx::Constraints *constr,
179 void update_sd_second_half(int64_t step,
180 real *dvdlambda, /* the contribution to be added to the bonded interactions */
181 const t_inputrec *inputrec, /* input record and box stuff */
186 gmx_wallcycle_t wcycle,
188 gmx::Constraints *constr,
192 void finish_update(const t_inputrec *inputrec,
195 const t_graph *graph,
197 gmx_wallcycle_t wcycle,
199 const gmx::Constraints *constr);
201 /* Return TRUE if OK, FALSE in case of Shake Error */
204 const rvec *x, const rvec *v, const matrix box,
205 const t_grpopts *opts, const t_mdatoms *md,
206 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel);
208 * Compute the partial kinetic energy for home particles;
209 * will be accumulated in the calling routine.
212 * Ekin = SUM(i) 0.5 m[i] v[i] (x) v[i]
214 * use v[i] = v[i] - u[i] when calculating temperature
216 * u must be accumulated already.
218 * Now also computes the contribution of the kinetic energy to the
225 init_ekinstate(ekinstate_t *ekinstate, const t_inputrec *ir);
228 update_ekinstate(ekinstate_t *ekinstate, const gmx_ekindata_t *ekind);
230 /*! \brief Restores data from \p ekinstate to \p ekind, then broadcasts it
231 to the rest of the simulation */
233 restore_ekinstate_from_state(const t_commrec *cr,
234 gmx_ekindata_t *ekind, const ekinstate_t *ekinstate);
236 void berendsen_tcoupl(const t_inputrec *ir, gmx_ekindata_t *ekind, real dt,
237 std::vector<double> &therm_integral); //NOLINT(google-runtime-references)
239 void andersen_tcoupl(const t_inputrec *ir, int64_t step,
240 const t_commrec *cr, const t_mdatoms *md,
241 gmx::ArrayRef<gmx::RVec> v,
242 real rate, const std::vector<bool> &randomize,
243 gmx::ArrayRef<const real> boltzfac);
245 void nosehoover_tcoupl(const t_grpopts *opts, const gmx_ekindata_t *ekind, real dt,
246 double xi[], double vxi[], const t_extmass *MassQ);
248 void trotter_update(const t_inputrec *ir, int64_t step, gmx_ekindata_t *ekind,
249 const gmx_enerdata_t *enerd, t_state *state, const tensor vir,
250 const t_mdatoms *md, const t_extmass *MassQ,
251 gmx::ArrayRef < std::vector < int>> trotter_seqlist, int trotter_seqno);
253 std::array < std::vector < int>, ettTSEQMAX> init_npt_vars(const t_inputrec *ir, t_state *state,
254 t_extmass *Mass, gmx_bool bTrotter);
256 real NPT_energy(const t_inputrec *ir, const t_state *state, const t_extmass *MassQ);
257 /* computes all the pressure/tempertature control energy terms to get a conserved energy */
259 void vrescale_tcoupl(const t_inputrec *ir, int64_t step,
260 gmx_ekindata_t *ekind, real dt,
261 double therm_integral[]);
262 /* Compute temperature scaling. For V-rescale it is done in update. */
264 void rescale_velocities(const gmx_ekindata_t *ekind, const t_mdatoms *mdatoms,
265 int start, int end, rvec v[]);
266 /* Rescale the velocities with the scaling factor in ekind */
268 //! Initialize simulated annealing.
269 bool initSimulatedAnnealing(t_inputrec *ir,
272 // TODO: This is the only function in update.h altering the inputrec
273 void update_annealing_target_temp(t_inputrec *ir, real t, gmx::Update *upd);
274 /* Set reference temp for simulated annealing at time t*/
276 real calc_temp(real ekin, real nrdf);
277 /* Calculate the temperature */
279 real calc_pres(int ePBC, int nwall, const matrix box, const tensor ekin, const tensor vir,
281 /* Calculate the pressure tensor, returns the scalar pressure.
282 * The unit of pressure is bar.
285 void parrinellorahman_pcoupl(FILE *fplog, int64_t step,
286 const t_inputrec *ir, real dt, const tensor pres,
287 const tensor box, tensor box_rel, tensor boxv,
289 gmx_bool bFirstStep);
291 void berendsen_pcoupl(FILE *fplog, int64_t step,
292 const t_inputrec *ir, real dt,
293 const tensor pres, const matrix box,
294 const matrix force_vir, const matrix constraint_vir,
295 matrix mu, double *baros_integral);
297 void berendsen_pscale(const t_inputrec *ir, const matrix mu,
298 matrix box, matrix box_rel,
299 int start, int nr_atoms,
300 rvec x[], const unsigned short cFREEZE[],
303 void pleaseCiteCouplingAlgorithms(FILE *fplog,
304 const t_inputrec &ir);
306 /*! \brief Computes the atom range for a thread to operate on, ensuring SIMD aligned ranges
308 * \param[in] numThreads The number of threads to divide atoms over
309 * \param[in] threadIndex The thread to get the range for
310 * \param[in] numAtoms The total number of atoms (on this rank)
311 * \param[out] startAtom The start of the atom range
312 * \param[out] endAtom The end of the atom range, note that this is in general not a multiple of the SIMD width
314 void getThreadAtomRange(int numThreads, int threadIndex, int numAtoms,
315 int *startAtom, int *endAtom);