Use std::vector for trotter_seq
[gromacs.git] / src / gromacs / mdlib / update.h
blobbe81265beac8f84f761c1718b1b4490afa16345d
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,2018, 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/real.h"
48 class ekinstate_t;
49 struct gmx_ekindata_t;
50 struct gmx_enerdata_t;
51 struct t_extmass;
52 struct t_fcdata;
53 struct t_graph;
54 struct t_grpopts;
55 struct t_inputrec;
56 struct t_mdatoms;
57 struct t_nrnb;
58 class t_state;
60 /* Abstract type for update */
61 struct gmx_update_t;
63 namespace gmx
65 class BoxDeformation;
66 class Constraints;
69 /* Initialize the stochastic dynamics struct */
70 gmx_update_t *init_update(const t_inputrec *ir,
71 gmx::BoxDeformation *deform);
73 /* Update pre-computed constants that depend on the reference
74 * temperature for coupling.
76 * This could change e.g. in simulated annealing. */
77 void update_temperature_constants(gmx_update_t *upd, const t_inputrec *ir);
79 /* Update the size of per-atom arrays (e.g. after DD re-partitioning,
80 which might increase the number of home atoms). */
81 void update_realloc(gmx_update_t *upd, int natoms);
83 /* Store the box at step step
84 * as a reference state for simulations with box deformation.
86 void set_deform_reference_box(gmx_update_t *upd,
87 int64_t step, matrix box);
89 void update_tcouple(int64_t step,
90 const t_inputrec *inputrec,
91 t_state *state,
92 gmx_ekindata_t *ekind,
93 const t_extmass *MassQ,
94 const t_mdatoms *md
97 /* Update Parrinello-Rahman, to be called before the coordinate update */
98 void update_pcouple_before_coordinates(FILE *fplog,
99 int64_t step,
100 const t_inputrec *inputrec,
101 t_state *state,
102 matrix parrinellorahmanMu,
103 matrix M,
104 gmx_bool bInitStep);
106 /* Update the box, to be called after the coordinate update.
107 * For Berendsen P-coupling, also calculates the scaling factor
108 * and scales the coordinates.
109 * When the deform option is used, scales coordinates and box here.
111 void update_pcouple_after_coordinates(FILE *fplog,
112 int64_t step,
113 const t_inputrec *inputrec,
114 const t_mdatoms *md,
115 const matrix pressure,
116 const matrix forceVirial,
117 const matrix constraintVirial,
118 const matrix parrinellorahmanMu,
119 t_state *state,
120 t_nrnb *nrnb,
121 gmx_update_t *upd);
123 void update_coords(int64_t step,
124 const t_inputrec *inputrec, /* input record and box stuff */
125 const t_mdatoms *md,
126 t_state *state,
127 gmx::ArrayRefWithPadding<gmx::RVec> f, /* forces on home particles */
128 const t_fcdata *fcd,
129 const gmx_ekindata_t *ekind,
130 const matrix M,
131 gmx_update_t *upd,
132 int bUpdatePart,
133 const t_commrec *cr, /* these shouldn't be here -- need to think about it */
134 const gmx::Constraints *constr);
136 /* Return TRUE if OK, FALSE in case of Shake Error */
138 extern gmx_bool update_randomize_velocities(const t_inputrec *ir, int64_t step, const t_commrec *cr,
139 const t_mdatoms *md,
140 gmx::ArrayRef<gmx::RVec> v,
141 const gmx_update_t *upd,
142 const gmx::Constraints *constr);
144 void constrain_velocities(int64_t step,
145 real *dvdlambda, /* the contribution to be added to the bonded interactions */
146 t_state *state,
147 tensor vir_part,
148 gmx::Constraints *constr,
149 gmx_bool bCalcVir,
150 bool do_log,
151 bool do_ene);
153 void constrain_coordinates(int64_t step,
154 real *dvdlambda, /* the contribution to be added to the bonded interactions */
155 t_state *state,
156 tensor vir_part,
157 gmx_update_t *upd,
158 gmx::Constraints *constr,
159 gmx_bool bCalcVir,
160 bool do_log,
161 bool do_ene);
163 void update_sd_second_half(int64_t step,
164 real *dvdlambda, /* the contribution to be added to the bonded interactions */
165 const t_inputrec *inputrec, /* input record and box stuff */
166 const t_mdatoms *md,
167 t_state *state,
168 const t_commrec *cr,
169 t_nrnb *nrnb,
170 gmx_wallcycle_t wcycle,
171 gmx_update_t *upd,
172 gmx::Constraints *constr,
173 bool do_log,
174 bool do_ene);
176 void finish_update(const t_inputrec *inputrec,
177 const t_mdatoms *md,
178 t_state *state,
179 const t_graph *graph,
180 t_nrnb *nrnb,
181 gmx_wallcycle_t wcycle,
182 gmx_update_t *upd,
183 const gmx::Constraints *constr);
185 /* Return TRUE if OK, FALSE in case of Shake Error */
187 void calc_ke_part(const t_state *state, const t_grpopts *opts, const t_mdatoms *md,
188 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel);
190 * Compute the partial kinetic energy for home particles;
191 * will be accumulated in the calling routine.
192 * The tensor is
194 * Ekin = SUM(i) 0.5 m[i] v[i] (x) v[i]
196 * use v[i] = v[i] - u[i] when calculating temperature
198 * u must be accumulated already.
200 * Now also computes the contribution of the kinetic energy to the
201 * free energy
206 void
207 init_ekinstate(ekinstate_t *ekinstate, const t_inputrec *ir);
209 void
210 update_ekinstate(ekinstate_t *ekinstate, const gmx_ekindata_t *ekind);
212 /*! \brief Restores data from \p ekinstate to \p ekind, then broadcasts it
213 to the rest of the simulation */
214 void
215 restore_ekinstate_from_state(const t_commrec *cr,
216 gmx_ekindata_t *ekind, const ekinstate_t *ekinstate);
218 void berendsen_tcoupl(const t_inputrec *ir, const gmx_ekindata_t *ekind, real dt,
219 std::vector<double> &therm_integral); //NOLINT(google-runtime-references)
221 void andersen_tcoupl(const t_inputrec *ir, int64_t step,
222 const t_commrec *cr, const t_mdatoms *md,
223 gmx::ArrayRef<gmx::RVec> v,
224 real rate, const gmx_bool *randomize, const real *boltzfac);
226 void nosehoover_tcoupl(const t_grpopts *opts, const gmx_ekindata_t *ekind, real dt,
227 double xi[], double vxi[], const t_extmass *MassQ);
229 void trotter_update(const t_inputrec *ir, int64_t step, gmx_ekindata_t *ekind,
230 const gmx_enerdata_t *enerd, t_state *state, const tensor vir,
231 const t_mdatoms *md, const t_extmass *MassQ,
232 gmx::ArrayRef < std::vector < int>> trotter_seqlist, int trotter_seqno);
234 std::array < std::vector < int>, ettTSEQMAX> init_npt_vars(const t_inputrec *ir, t_state *state,
235 t_extmass *Mass, gmx_bool bTrotter);
237 real NPT_energy(const t_inputrec *ir, const t_state *state, const t_extmass *MassQ);
238 /* computes all the pressure/tempertature control energy terms to get a conserved energy */
240 void vrescale_tcoupl(const t_inputrec *ir, int64_t step,
241 gmx_ekindata_t *ekind, real dt,
242 double therm_integral[]);
243 /* Compute temperature scaling. For V-rescale it is done in update. */
245 void rescale_velocities(const gmx_ekindata_t *ekind, const t_mdatoms *mdatoms,
246 int start, int end, rvec v[]);
247 /* Rescale the velocities with the scaling factor in ekind */
249 // TODO: This is the only function in update.h altering the inputrec
250 void update_annealing_target_temp(t_inputrec *ir, real t, gmx_update_t *upd);
251 /* Set reference temp for simulated annealing at time t*/
253 real calc_temp(real ekin, real nrdf);
254 /* Calculate the temperature */
256 real calc_pres(int ePBC, int nwall, const matrix box, const tensor ekin, const tensor vir,
257 tensor pres);
258 /* Calculate the pressure tensor, returns the scalar pressure.
259 * The unit of pressure is bar.
262 void parrinellorahman_pcoupl(FILE *fplog, int64_t step,
263 const t_inputrec *ir, real dt, const tensor pres,
264 const tensor box, tensor box_rel, tensor boxv,
265 tensor M, matrix mu,
266 gmx_bool bFirstStep);
268 void berendsen_pcoupl(FILE *fplog, int64_t step,
269 const t_inputrec *ir, real dt,
270 const tensor pres, const matrix box,
271 const matrix force_vir, const matrix constraint_vir,
272 matrix mu, double *baros_integral);
274 void berendsen_pscale(const t_inputrec *ir, const matrix mu,
275 matrix box, matrix box_rel,
276 int start, int nr_atoms,
277 rvec x[], const unsigned short cFREEZE[],
278 t_nrnb *nrnb);
279 #endif