Merge branch 'release-2020' into master
[alexxy/gromacs.git] / src / gromacs / mdlib / coupling.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 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.
11  *
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.
16  *
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.
21  *
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.
26  *
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.
34  *
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.
37  */
38 #ifndef GMX_MDLIB_COUPLING_H
39 #define GMX_MDLIB_COUPLING_H
40
41 #include "gromacs/math/paddedvector.h"
42 #include "gromacs/math/vectypes.h"
43 #include "gromacs/mdtypes/md_enums.h"
44 #include "gromacs/utility/arrayref.h"
45 #include "gromacs/utility/basedefinitions.h"
46 #include "gromacs/utility/real.h"
47
48 struct gmx_ekindata_t;
49 struct gmx_enerdata_t;
50 struct t_commrec;
51 struct t_extmass;
52 struct t_grpopts;
53 struct t_inputrec;
54 struct t_mdatoms;
55 struct t_nrnb;
56 class t_state;
57
58 enum class PbcType;
59
60 namespace gmx
61 {
62 class BoxDeformation;
63 class Constraints;
64 class Update;
65 }; // namespace gmx
66
67 /* Update the size of per-atom arrays (e.g. after DD re-partitioning,
68    which might increase the number of home atoms). */
69
70 void update_tcouple(int64_t           step,
71                     const t_inputrec* inputrec,
72                     t_state*          state,
73                     gmx_ekindata_t*   ekind,
74                     const t_extmass*  MassQ,
75                     const t_mdatoms*  md);
76
77 /* Update Parrinello-Rahman, to be called before the coordinate update */
78 void update_pcouple_before_coordinates(FILE*             fplog,
79                                        int64_t           step,
80                                        const t_inputrec* inputrec,
81                                        t_state*          state,
82                                        matrix            parrinellorahmanMu,
83                                        matrix            M,
84                                        gmx_bool          bInitStep);
85
86 /* Update the box, to be called after the coordinate update.
87  * For Berendsen P-coupling, also calculates the scaling factor
88  * and scales the coordinates.
89  * When the deform option is used, scales coordinates and box here.
90  */
91 void update_pcouple_after_coordinates(FILE*                fplog,
92                                       int64_t              step,
93                                       const t_inputrec*    inputrec,
94                                       const t_mdatoms*     md,
95                                       const matrix         pressure,
96                                       const matrix         forceVirial,
97                                       const matrix         constraintVirial,
98                                       matrix               pressureCouplingMu,
99                                       t_state*             state,
100                                       t_nrnb*              nrnb,
101                                       gmx::BoxDeformation* boxDeformation,
102                                       bool                 scaleCoordinates);
103
104 /* Return TRUE if OK, FALSE in case of Shake Error */
105
106 extern gmx_bool update_randomize_velocities(const t_inputrec*        ir,
107                                             int64_t                  step,
108                                             const t_commrec*         cr,
109                                             const t_mdatoms*         md,
110                                             gmx::ArrayRef<gmx::RVec> v,
111                                             const gmx::Update*       upd,
112                                             const gmx::Constraints*  constr);
113
114 void berendsen_tcoupl(const t_inputrec*    ir,
115                       gmx_ekindata_t*      ekind,
116                       real                 dt,
117                       std::vector<double>& therm_integral); //NOLINT(google-runtime-references)
118
119 void andersen_tcoupl(const t_inputrec*         ir,
120                      int64_t                   step,
121                      const t_commrec*          cr,
122                      const t_mdatoms*          md,
123                      gmx::ArrayRef<gmx::RVec>  v,
124                      real                      rate,
125                      const std::vector<bool>&  randomize,
126                      gmx::ArrayRef<const real> boltzfac);
127
128 void nosehoover_tcoupl(const t_grpopts*      opts,
129                        const gmx_ekindata_t* ekind,
130                        real                  dt,
131                        double                xi[],
132                        double                vxi[],
133                        const t_extmass*      MassQ);
134
135 void trotter_update(const t_inputrec*               ir,
136                     int64_t                         step,
137                     gmx_ekindata_t*                 ekind,
138                     const gmx_enerdata_t*           enerd,
139                     t_state*                        state,
140                     const tensor                    vir,
141                     const t_mdatoms*                md,
142                     const t_extmass*                MassQ,
143                     gmx::ArrayRef<std::vector<int>> trotter_seqlist,
144                     int                             trotter_seqno);
145
146 std::array<std::vector<int>, ettTSEQMAX>
147 init_npt_vars(const t_inputrec* ir, t_state* state, t_extmass* Mass, gmx_bool bTrotter);
148
149 real NPT_energy(const t_inputrec* ir, const t_state* state, const t_extmass* MassQ);
150 /* computes all the pressure/tempertature control energy terms to get a conserved energy */
151
152 void vrescale_tcoupl(const t_inputrec* ir, int64_t step, gmx_ekindata_t* ekind, real dt, double therm_integral[]);
153 /* Compute temperature scaling. For V-rescale it is done in update. */
154
155 void rescale_velocities(const gmx_ekindata_t* ekind, const t_mdatoms* mdatoms, int start, int end, rvec v[]);
156 /* Rescale the velocities with the scaling factor in ekind */
157
158 //! Check whether we do simulated annealing.
159 bool doSimulatedAnnealing(const t_inputrec* ir);
160
161 //! Initialize simulated annealing.
162 bool initSimulatedAnnealing(t_inputrec* ir, gmx::Update* upd);
163
164 // TODO: This is the only function in update.h altering the inputrec
165 void update_annealing_target_temp(t_inputrec* ir, real t, gmx::Update* upd);
166 /* Set reference temp for simulated annealing at time t*/
167
168 real calc_temp(real ekin, real nrdf);
169 /* Calculate the temperature */
170
171 real calc_pres(PbcType pbcType, int nwall, const matrix box, const tensor ekin, const tensor vir, tensor pres);
172 /* Calculate the pressure tensor, returns the scalar pressure.
173  * The unit of pressure is bar.
174  */
175
176 void parrinellorahman_pcoupl(FILE*             fplog,
177                              int64_t           step,
178                              const t_inputrec* ir,
179                              real              dt,
180                              const tensor      pres,
181                              const tensor      box,
182                              tensor            box_rel,
183                              tensor            boxv,
184                              tensor            M,
185                              matrix            mu,
186                              gmx_bool          bFirstStep);
187
188 void berendsen_pcoupl(FILE*             fplog,
189                       int64_t           step,
190                       const t_inputrec* ir,
191                       real              dt,
192                       const tensor      pres,
193                       const matrix      box,
194                       const matrix      force_vir,
195                       const matrix      constraint_vir,
196                       matrix            mu,
197                       double*           baros_integral);
198
199 void berendsen_pscale(const t_inputrec*    ir,
200                       const matrix         mu,
201                       matrix               box,
202                       matrix               box_rel,
203                       int                  start,
204                       int                  nr_atoms,
205                       rvec                 x[],
206                       const unsigned short cFREEZE[],
207                       t_nrnb*              nrnb,
208                       bool                 scaleCoordinates);
209
210 void pleaseCiteCouplingAlgorithms(FILE* fplog, const t_inputrec& ir);
211
212 /*! \brief Generate a new kinetic energy for the v-rescale thermostat
213  *
214  * Generates a new value for the kinetic energy, according to
215  * Bussi et al JCP (2007), Eq. (A7)
216  *
217  * This is used by update_tcoupl(), and by the VRescaleThermostat of the modular
218  * simulator.
219  * \todo Move this to the VRescaleThermostat once the modular simulator becomes
220  *       the default code path.
221  *
222  * \param[in] kk     present value of the kinetic energy of the atoms to be thermalized (in
223  * arbitrary units) \param[in] sigma  target average value of the kinetic energy (ndeg k_b T/2)  (in
224  * the same units as kk) \param[in] ndeg   number of degrees of freedom of the atoms to be
225  * thermalized \param[in] taut   relaxation time of the thermostat, in units of 'how often this
226  * routine is called' \param[in] step   the time step this routine is called on \param[in] seed the
227  * random number generator seed \return  the new kinetic energy
228  */
229 real vrescale_resamplekin(real kk, real sigma, real ndeg, real taut, int64_t step, int64_t seed);
230
231 #endif // GMX_MDLIB_COUPLING_H