Apply clang-format to source tree
[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                                       const matrix      parrinellorahmanMu,
135                                       t_state*          state,
136                                       t_nrnb*           nrnb,
137                                       gmx::Update*      upd);
138
139 void update_coords(int64_t           step,
140                    const t_inputrec* inputrec, /* input record and box stuff    */
141                    const t_mdatoms*  md,
142                    t_state*          state,
143                    gmx::ArrayRefWithPadding<const gmx::RVec> f, /* forces on home particles */
144                    const t_fcdata*                           fcd,
145                    const gmx_ekindata_t*                     ekind,
146                    const matrix                              M,
147                    gmx::Update*                              upd,
148                    int                                       bUpdatePart,
149                    const t_commrec* cr, /* these shouldn't be here -- need to think about it */
150                    const gmx::Constraints* constr);
151
152 /* Return TRUE if OK, FALSE in case of Shake Error */
153
154 extern gmx_bool update_randomize_velocities(const t_inputrec*        ir,
155                                             int64_t                  step,
156                                             const t_commrec*         cr,
157                                             const t_mdatoms*         md,
158                                             gmx::ArrayRef<gmx::RVec> v,
159                                             const gmx::Update*       upd,
160                                             const gmx::Constraints*  constr);
161
162 void constrain_velocities(int64_t step,
163                           real* dvdlambda, /* the contribution to be added to the bonded interactions */
164                           t_state*          state,
165                           tensor            vir_part,
166                           gmx::Constraints* constr,
167                           gmx_bool          bCalcVir,
168                           bool              do_log,
169                           bool              do_ene);
170
171 void constrain_coordinates(int64_t step,
172                            real* dvdlambda, /* the contribution to be added to the bonded interactions */
173                            t_state*          state,
174                            tensor            vir_part,
175                            gmx::Update*      upd,
176                            gmx::Constraints* constr,
177                            gmx_bool          bCalcVir,
178                            bool              do_log,
179                            bool              do_ene);
180
181 void update_sd_second_half(int64_t step,
182                            real* dvdlambda, /* the contribution to be added to the bonded interactions */
183                            const t_inputrec* inputrec, /* input record and box stuff */
184                            const t_mdatoms*  md,
185                            t_state*          state,
186                            const t_commrec*  cr,
187                            t_nrnb*           nrnb,
188                            gmx_wallcycle_t   wcycle,
189                            gmx::Update*      upd,
190                            gmx::Constraints* constr,
191                            bool              do_log,
192                            bool              do_ene);
193
194 void finish_update(const t_inputrec*       inputrec,
195                    const t_mdatoms*        md,
196                    t_state*                state,
197                    const t_graph*          graph,
198                    t_nrnb*                 nrnb,
199                    gmx_wallcycle_t         wcycle,
200                    gmx::Update*            upd,
201                    const gmx::Constraints* constr);
202
203 /* Return TRUE if OK, FALSE in case of Shake Error */
204
205 void calc_ke_part(const rvec*      x,
206                   const rvec*      v,
207                   const matrix     box,
208                   const t_grpopts* opts,
209                   const t_mdatoms* md,
210                   gmx_ekindata_t*  ekind,
211                   t_nrnb*          nrnb,
212                   gmx_bool         bEkinAveVel);
213 /*
214  * Compute the partial kinetic energy for home particles;
215  * will be accumulated in the calling routine.
216  * The tensor is
217  *
218  * Ekin = SUM(i) 0.5 m[i] v[i] (x) v[i]
219  *
220  *     use v[i] = v[i] - u[i] when calculating temperature
221  *
222  * u must be accumulated already.
223  *
224  * Now also computes the contribution of the kinetic energy to the
225  * free energy
226  *
227  */
228
229
230 void init_ekinstate(ekinstate_t* ekinstate, const t_inputrec* ir);
231
232 void update_ekinstate(ekinstate_t* ekinstate, const gmx_ekindata_t* ekind);
233
234 /*! \brief Restores data from \p ekinstate to \p ekind, then broadcasts it
235    to the rest of the simulation */
236 void restore_ekinstate_from_state(const t_commrec* cr, gmx_ekindata_t* ekind, const ekinstate_t* ekinstate);
237
238 void berendsen_tcoupl(const t_inputrec*    ir,
239                       gmx_ekindata_t*      ekind,
240                       real                 dt,
241                       std::vector<double>& therm_integral); //NOLINT(google-runtime-references)
242
243 void andersen_tcoupl(const t_inputrec*         ir,
244                      int64_t                   step,
245                      const t_commrec*          cr,
246                      const t_mdatoms*          md,
247                      gmx::ArrayRef<gmx::RVec>  v,
248                      real                      rate,
249                      const std::vector<bool>&  randomize,
250                      gmx::ArrayRef<const real> boltzfac);
251
252 void nosehoover_tcoupl(const t_grpopts*      opts,
253                        const gmx_ekindata_t* ekind,
254                        real                  dt,
255                        double                xi[],
256                        double                vxi[],
257                        const t_extmass*      MassQ);
258
259 void trotter_update(const t_inputrec*               ir,
260                     int64_t                         step,
261                     gmx_ekindata_t*                 ekind,
262                     const gmx_enerdata_t*           enerd,
263                     t_state*                        state,
264                     const tensor                    vir,
265                     const t_mdatoms*                md,
266                     const t_extmass*                MassQ,
267                     gmx::ArrayRef<std::vector<int>> trotter_seqlist,
268                     int                             trotter_seqno);
269
270 std::array<std::vector<int>, ettTSEQMAX>
271 init_npt_vars(const t_inputrec* ir, t_state* state, t_extmass* Mass, gmx_bool bTrotter);
272
273 real NPT_energy(const t_inputrec* ir, const t_state* state, const t_extmass* MassQ);
274 /* computes all the pressure/tempertature control energy terms to get a conserved energy */
275
276 void vrescale_tcoupl(const t_inputrec* ir, int64_t step, gmx_ekindata_t* ekind, real dt, double therm_integral[]);
277 /* Compute temperature scaling. For V-rescale it is done in update. */
278
279 void rescale_velocities(const gmx_ekindata_t* ekind, const t_mdatoms* mdatoms, int start, int end, rvec v[]);
280 /* Rescale the velocities with the scaling factor in ekind */
281
282 //! Check whether we do simulated annealing.
283 bool doSimulatedAnnealing(const t_inputrec* ir);
284
285 //! Initialize simulated annealing.
286 bool initSimulatedAnnealing(t_inputrec* ir, gmx::Update* upd);
287
288 // TODO: This is the only function in update.h altering the inputrec
289 void update_annealing_target_temp(t_inputrec* ir, real t, gmx::Update* upd);
290 /* Set reference temp for simulated annealing at time t*/
291
292 real calc_temp(real ekin, real nrdf);
293 /* Calculate the temperature */
294
295 real calc_pres(int ePBC, int nwall, const matrix box, const tensor ekin, const tensor vir, tensor pres);
296 /* Calculate the pressure tensor, returns the scalar pressure.
297  * The unit of pressure is bar.
298  */
299
300 void parrinellorahman_pcoupl(FILE*             fplog,
301                              int64_t           step,
302                              const t_inputrec* ir,
303                              real              dt,
304                              const tensor      pres,
305                              const tensor      box,
306                              tensor            box_rel,
307                              tensor            boxv,
308                              tensor            M,
309                              matrix            mu,
310                              gmx_bool          bFirstStep);
311
312 void berendsen_pcoupl(FILE*             fplog,
313                       int64_t           step,
314                       const t_inputrec* ir,
315                       real              dt,
316                       const tensor      pres,
317                       const matrix      box,
318                       const matrix      force_vir,
319                       const matrix      constraint_vir,
320                       matrix            mu,
321                       double*           baros_integral);
322
323 void berendsen_pscale(const t_inputrec*    ir,
324                       const matrix         mu,
325                       matrix               box,
326                       matrix               box_rel,
327                       int                  start,
328                       int                  nr_atoms,
329                       rvec                 x[],
330                       const unsigned short cFREEZE[],
331                       t_nrnb*              nrnb);
332
333 void pleaseCiteCouplingAlgorithms(FILE* fplog, const t_inputrec& ir);
334
335 /*! \brief Computes the atom range for a thread to operate on, ensuring SIMD aligned ranges
336  *
337  * \param[in]  numThreads   The number of threads to divide atoms over
338  * \param[in]  threadIndex  The thread to get the range for
339  * \param[in]  numAtoms     The total number of atoms (on this rank)
340  * \param[out] startAtom    The start of the atom range
341  * \param[out] endAtom      The end of the atom range, note that this is in general not a multiple of the SIMD width
342  */
343 void getThreadAtomRange(int numThreads, int threadIndex, int numAtoms, int* startAtom, int* endAtom);
344
345 /*! \brief Generate a new kinetic energy for the v-rescale thermostat
346  *
347  * Generates a new value for the kinetic energy, according to
348  * Bussi et al JCP (2007), Eq. (A7)
349  *
350  * This is used by update_tcoupl(), and by the VRescaleThermostat of the modular
351  * simulator.
352  * TODO: Move this to the VRescaleThermostat once the modular simulator becomes
353  *       the default code path.
354  *
355  * @param kk     present value of the kinetic energy of the atoms to be thermalized (in arbitrary units)
356  * @param sigma  target average value of the kinetic energy (ndeg k_b T/2)  (in the same units as kk)
357  * @param ndeg   number of degrees of freedom of the atoms to be thermalized
358  * @param taut   relaxation time of the thermostat, in units of 'how often this routine is called'
359  * @param step   the time step this routine is called on
360  * @param seed   the random number generator seed
361  * @return  the new kinetic energy
362  */
363 real vrescale_resamplekin(real kk, real sigma, real ndeg, real taut, int64_t step, int64_t seed);
364
365 #endif