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