Make PBC type enumeration into PbcType enum class
[gromacs.git] / src / gromacs / mdlib / update.h
blob5b9607b2aa37a18122fc5bac6b122c960287ea2e
1 /*
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"
50 class ekinstate_t;
51 struct gmx_ekindata_t;
52 struct gmx_enerdata_t;
53 enum class PbcType;
54 struct t_extmass;
55 struct t_fcdata;
56 struct t_graph;
57 struct t_grpopts;
58 struct t_inputrec;
59 struct t_mdatoms;
60 struct t_nrnb;
61 class t_state;
63 /* Abstract type for update */
64 struct gmx_stochd_t;
66 namespace gmx
68 class BoxDeformation;
69 class Constraints;
72 /*! \libinternal
73 * \brief Contains data for update phase */
74 class Update
76 public:
77 //! Constructor
78 Update(const t_inputrec* ir, BoxDeformation* boxDeformation);
79 ~Update();
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;
87 //! Resizes xp
88 void setNumAtoms(int nAtoms);
90 private:
91 //! Implementation type.
92 class Impl;
93 //! Implementation object.
94 PrivateImplPointer<Impl> impl_;
97 }; // namespace gmx
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,
110 t_state* state,
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,
117 int64_t step,
118 const t_inputrec* inputrec,
119 t_state* state,
120 matrix parrinellorahmanMu,
121 matrix M,
122 gmx_bool bInitStep);
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,
130 int64_t step,
131 const t_inputrec* inputrec,
132 const t_mdatoms* md,
133 const matrix pressure,
134 const matrix forceVirial,
135 const matrix constraintVirial,
136 matrix pressureCouplingMu,
137 t_state* state,
138 t_nrnb* nrnb,
139 gmx::Update* upd,
140 bool scaleCoordinates);
142 void update_coords(int64_t step,
143 const t_inputrec* inputrec, /* input record and box stuff */
144 const t_mdatoms* md,
145 t_state* state,
146 gmx::ArrayRefWithPadding<const gmx::RVec> f, /* forces on home particles */
147 const t_fcdata* fcd,
148 const gmx_ekindata_t* ekind,
149 const matrix M,
150 gmx::Update* upd,
151 int bUpdatePart,
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,
158 int64_t step,
159 const t_commrec* cr,
160 const t_mdatoms* md,
161 gmx::ArrayRef<gmx::RVec> v,
162 const gmx::Update* upd,
163 const gmx::Constraints* constr);
165 void constrain_velocities(int64_t step,
166 real* dvdlambda, /* the contribution to be added to the bonded interactions */
167 t_state* state,
168 tensor vir_part,
169 gmx::Constraints* constr,
170 gmx_bool bCalcVir,
171 bool do_log,
172 bool do_ene);
174 void constrain_coordinates(int64_t step,
175 real* dvdlambda, /* the contribution to be added to the bonded interactions */
176 t_state* state,
177 tensor vir_part,
178 gmx::Update* upd,
179 gmx::Constraints* constr,
180 gmx_bool bCalcVir,
181 bool do_log,
182 bool do_ene);
184 void update_sd_second_half(int64_t step,
185 real* dvdlambda, /* the contribution to be added to the bonded interactions */
186 const t_inputrec* inputrec, /* input record and box stuff */
187 const t_mdatoms* md,
188 t_state* state,
189 const t_commrec* cr,
190 t_nrnb* nrnb,
191 gmx_wallcycle_t wcycle,
192 gmx::Update* upd,
193 gmx::Constraints* constr,
194 bool do_log,
195 bool do_ene);
197 void finish_update(const t_inputrec* inputrec,
198 const t_mdatoms* md,
199 t_state* state,
200 const t_graph* graph,
201 t_nrnb* nrnb,
202 gmx_wallcycle_t wcycle,
203 gmx::Update* upd,
204 const gmx::Constraints* constr);
206 /* Return TRUE if OK, FALSE in case of Shake Error */
208 void calc_ke_part(const rvec* x,
209 const rvec* v,
210 const matrix box,
211 const t_grpopts* opts,
212 const t_mdatoms* md,
213 gmx_ekindata_t* ekind,
214 t_nrnb* nrnb,
215 gmx_bool bEkinAveVel);
217 * Compute the partial kinetic energy for home particles;
218 * will be accumulated in the calling routine.
219 * The tensor is
221 * Ekin = SUM(i) 0.5 m[i] v[i] (x) v[i]
223 * use v[i] = v[i] - u[i] when calculating temperature
225 * u must be accumulated already.
227 * Now also computes the contribution of the kinetic energy to the
228 * free energy
233 void init_ekinstate(ekinstate_t* ekinstate, const t_inputrec* ir);
235 void update_ekinstate(ekinstate_t* ekinstate, const gmx_ekindata_t* ekind);
237 /*! \brief Restores data from \p ekinstate to \p ekind, then broadcasts it
238 to the rest of the simulation */
239 void restore_ekinstate_from_state(const t_commrec* cr, gmx_ekindata_t* ekind, const ekinstate_t* ekinstate);
241 void berendsen_tcoupl(const t_inputrec* ir,
242 gmx_ekindata_t* ekind,
243 real dt,
244 std::vector<double>& therm_integral); //NOLINT(google-runtime-references)
246 void andersen_tcoupl(const t_inputrec* ir,
247 int64_t step,
248 const t_commrec* cr,
249 const t_mdatoms* md,
250 gmx::ArrayRef<gmx::RVec> v,
251 real rate,
252 const std::vector<bool>& randomize,
253 gmx::ArrayRef<const real> boltzfac);
255 void nosehoover_tcoupl(const t_grpopts* opts,
256 const gmx_ekindata_t* ekind,
257 real dt,
258 double xi[],
259 double vxi[],
260 const t_extmass* MassQ);
262 void trotter_update(const t_inputrec* ir,
263 int64_t step,
264 gmx_ekindata_t* ekind,
265 const gmx_enerdata_t* enerd,
266 t_state* state,
267 const tensor vir,
268 const t_mdatoms* md,
269 const t_extmass* MassQ,
270 gmx::ArrayRef<std::vector<int>> trotter_seqlist,
271 int trotter_seqno);
273 std::array<std::vector<int>, ettTSEQMAX>
274 init_npt_vars(const t_inputrec* ir, t_state* state, t_extmass* Mass, gmx_bool bTrotter);
276 real NPT_energy(const t_inputrec* ir, const t_state* state, const t_extmass* MassQ);
277 /* computes all the pressure/tempertature control energy terms to get a conserved energy */
279 void vrescale_tcoupl(const t_inputrec* ir, int64_t step, gmx_ekindata_t* ekind, real dt, double therm_integral[]);
280 /* Compute temperature scaling. For V-rescale it is done in update. */
282 void rescale_velocities(const gmx_ekindata_t* ekind, const t_mdatoms* mdatoms, int start, int end, rvec v[]);
283 /* Rescale the velocities with the scaling factor in ekind */
285 //! Check whether we do simulated annealing.
286 bool doSimulatedAnnealing(const t_inputrec* ir);
288 //! Initialize simulated annealing.
289 bool initSimulatedAnnealing(t_inputrec* ir, gmx::Update* upd);
291 // TODO: This is the only function in update.h altering the inputrec
292 void update_annealing_target_temp(t_inputrec* ir, real t, gmx::Update* upd);
293 /* Set reference temp for simulated annealing at time t*/
295 real calc_temp(real ekin, real nrdf);
296 /* Calculate the temperature */
298 real calc_pres(PbcType pbcType, int nwall, const matrix box, const tensor ekin, const tensor vir, tensor pres);
299 /* Calculate the pressure tensor, returns the scalar pressure.
300 * The unit of pressure is bar.
303 void parrinellorahman_pcoupl(FILE* fplog,
304 int64_t step,
305 const t_inputrec* ir,
306 real dt,
307 const tensor pres,
308 const tensor box,
309 tensor box_rel,
310 tensor boxv,
311 tensor M,
312 matrix mu,
313 gmx_bool bFirstStep);
315 void berendsen_pcoupl(FILE* fplog,
316 int64_t step,
317 const t_inputrec* ir,
318 real dt,
319 const tensor pres,
320 const matrix box,
321 const matrix force_vir,
322 const matrix constraint_vir,
323 matrix mu,
324 double* baros_integral);
326 void berendsen_pscale(const t_inputrec* ir,
327 const matrix mu,
328 matrix box,
329 matrix box_rel,
330 int start,
331 int nr_atoms,
332 rvec x[],
333 const unsigned short cFREEZE[],
334 t_nrnb* nrnb,
335 bool scaleCoordinates);
337 void pleaseCiteCouplingAlgorithms(FILE* fplog, const t_inputrec& ir);
339 /*! \brief Computes the atom range for a thread to operate on, ensuring SIMD aligned ranges
341 * \param[in] numThreads The number of threads to divide atoms over
342 * \param[in] threadIndex The thread to get the range for
343 * \param[in] numAtoms The total number of atoms (on this rank)
344 * \param[out] startAtom The start of the atom range
345 * \param[out] endAtom The end of the atom range, note that this is in general not a multiple of the SIMD width
347 void getThreadAtomRange(int numThreads, int threadIndex, int numAtoms, int* startAtom, int* endAtom);
349 /*! \brief Generate a new kinetic energy for the v-rescale thermostat
351 * Generates a new value for the kinetic energy, according to
352 * Bussi et al JCP (2007), Eq. (A7)
354 * This is used by update_tcoupl(), and by the VRescaleThermostat of the modular
355 * simulator.
356 * TODO: Move this to the VRescaleThermostat once the modular simulator becomes
357 * the default code path.
359 * @param kk present value of the kinetic energy of the atoms to be thermalized (in arbitrary units)
360 * @param sigma target average value of the kinetic energy (ndeg k_b T/2) (in the same units as kk)
361 * @param ndeg number of degrees of freedom of the atoms to be thermalized
362 * @param taut relaxation time of the thermostat, in units of 'how often this routine is called'
363 * @param step the time step this routine is called on
364 * @param seed the random number generator seed
365 * @return the new kinetic energy
367 real vrescale_resamplekin(real kk, real sigma, real ndeg, real taut, int64_t step, int64_t seed);
369 #endif