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