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