Constructor for gmx_ekindata_t
[alexxy/gromacs.git] / src / gromacs / mdlib / update_vv.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_UPDATE_VV_H
39 #define GMX_MDLIB_UPDATE_VV_H
40
41 #include <vector>
42
43 #include "gromacs/math/vectypes.h"
44 #include "gromacs/mdrunutility/handlerestart.h"
45 #include "gromacs/mdtypes/md_enums.h"
46
47 class gmx_ekindata_t;
48 struct gmx_enerdata_t;
49 struct gmx_global_stat;
50 struct gmx_localtop_t;
51 struct gmx_mtop_t;
52 struct gmx_wallcycle;
53 struct pull_t;
54 struct t_commrec;
55 struct t_extmass;
56 struct t_fcdata;
57 struct t_forcerec;
58 struct t_inputrec;
59 struct t_mdatoms;
60 struct t_nrnb;
61 class t_state;
62 struct t_vcm;
63
64 namespace gmx
65 {
66 class Constraints;
67 class ForceBuffers;
68 class MDLogger;
69 class SimulationSignaller;
70 class Update;
71 } // namespace gmx
72
73 /*! \brief Make the first step of Velocity Verlet integration
74  *
75  * \param[in]  step              Current timestep.
76  * \param[in]  bFirstStep        Is it a first step.
77  * \param[in]  bInitStep         Is it an initialization step.
78  * \param[in]  startingBehavior  Describes whether this is a restart appending to output files.
79  * \param[in]  nstglobalcomm     Will globals be computed on this step.
80  * \param[in]  ir                Input record.
81  * \param[in]  fr                Force record.
82  * \param[in]  cr                Comunication record.
83  * \param[in]  state             Simulation state.
84  * \param[in]  mdatoms           MD atoms data.
85  * \param[in]  fcdata            Force calculation data.
86  * \param[in]  MassQ             Mass/pressure data.
87  * \param[in]  vcm               Center of mass motion removal.
88  * \param[in]  top_global        Global topology.
89  * \param[in]  top               Local topology.
90  * \param[in]  enerd             Energy data.
91  * \param[in]  ekind             Kinetic energy data.
92  * \param[in]  gstat             Storage of thermodynamic parameters data.
93  * \param[out] last_ekin         Kinetic energies of the last step.
94  * \param[in]  bCalcVir          If the virial is computed on this step.
95  * \param[in]  total_vir         Total virial tensor.
96  * \param[in]  shake_vir         Constraints virial.
97  * \param[in]  force_vir         Force virial.
98  * \param[in]  pres              Pressure tensor.
99  * \param[in]  M                 Parrinello-Rahman velocity scaling matrix.
100  * \param[in]  do_log            Do logging on this step.
101  * \param[in]  do_ene            Print energies on this step.
102  * \param[in]  bCalcEner         Compute energies on this step.
103  * \param[in]  bGStat            Collect globals this step.
104  * \param[in]  bStopCM           Stop the center of mass motion on this step.
105  * \param[in]  bTrotter          Do trotter routines this step.
106  * \param[in]  bExchanged        If this is a replica exchange step.
107  * \param[out] bSumEkinhOld      Old kinetic energies will need to be summed up.
108  * \param[out] saved_conserved_quantity  Place to store the conserved energy.
109  * \param[in]  f                 Force buffers.
110  * \param[in]  upd               Update object.
111  * \param[in]  constr            Constraints object.
112  * \param[in]  nullSignaller     Simulation signaller.
113  * \param[in]  trotter_seq       NPT variables.
114  * \param[in]  nrnb              Cycle counters.
115  * \param[in]  mdlog             Logger.
116  * \param[in]  fplog             Another logger.
117  * \param[in]  wcycle            Wall-clock cycle counter.
118  */
119 void integrateVVFirstStep(int64_t                                  step,
120                           bool                                     bFirstStep,
121                           bool                                     bInitStep,
122                           gmx::StartingBehavior                    startingBehavior,
123                           int                                      nstglobalcomm,
124                           const t_inputrec*                        ir,
125                           t_forcerec*                              fr,
126                           t_commrec*                               cr,
127                           t_state*                                 state,
128                           t_mdatoms*                               mdatoms,
129                           const t_fcdata&                          fcdata,
130                           t_extmass*                               MassQ,
131                           t_vcm*                                   vcm,
132                           const gmx_mtop_t&                        top_global,
133                           const gmx_localtop_t&                    top,
134                           gmx_enerdata_t*                          enerd,
135                           gmx_ekindata_t*                          ekind,
136                           gmx_global_stat*                         gstat,
137                           real*                                    last_ekin,
138                           bool                                     bCalcVir,
139                           tensor                                   total_vir,
140                           tensor                                   shake_vir,
141                           tensor                                   force_vir,
142                           tensor                                   pres,
143                           matrix                                   M,
144                           bool                                     do_log,
145                           bool                                     do_ene,
146                           bool                                     bCalcEner,
147                           bool                                     bGStat,
148                           bool                                     bStopCM,
149                           bool                                     bTrotter,
150                           bool                                     bExchanged,
151                           bool*                                    bSumEkinhOld,
152                           real*                                    saved_conserved_quantity,
153                           gmx::ForceBuffers*                       f,
154                           gmx::Update*                             upd,
155                           gmx::Constraints*                        constr,
156                           gmx::SimulationSignaller*                nullSignaller,
157                           std::array<std::vector<int>, ettTSEQMAX> trotter_seq,
158                           t_nrnb*                                  nrnb,
159                           const gmx::MDLogger&                     mdlog,
160                           FILE*                                    fplog,
161                           gmx_wallcycle*                           wcycle);
162
163
164 /*! \brief Make the second step of Velocity Verlet integration
165  *
166  * \param[in]  step              Current timestep.
167  * \param[in]  ir                Input record.
168  * \param[in]  fr                Force record.
169  * \param[in]  cr                Comunication record.
170  * \param[in]  state             Simulation state.
171  * \param[in]  mdatoms           MD atoms data.
172  * \param[in]  fcdata            Force calculation data.
173  * \param[in]  MassQ             Mass/pressure data.
174  * \param[in]  vcm               Center of mass motion removal.
175  * \param[in]  pull_work         Pulling data.
176  * \param[in]  enerd             Energy data.
177  * \param[in]  ekind             Kinetic energy data.
178  * \param[in]  gstat             Storage of thermodynamic parameters data.
179  * \param[out] dvdl_constr       FEP data for constraints.
180  * \param[in]  bCalcVir          If the virial is computed on this step.
181  * \param[in]  total_vir         Total virial tensor.
182  * \param[in]  shake_vir         Constraints virial.
183  * \param[in]  force_vir         Force virial.
184  * \param[in]  pres              Pressure tensor.
185  * \param[in]  M                 Parrinello-Rahman velocity scaling matrix.
186  * \param[in]  lastbox           Last recorded PBC box.
187  * \param[in]  do_log            Do logging on this step.
188  * \param[in]  do_ene            Print energies on this step.
189  * \param[in]  bGStat            Collect globals this step.
190  * \param[out] bSumEkinhOld      Old kinetic energies need to be summed up.
191  * \param[in]  f                 Force buffers.
192  * \param[in]  cbuf              Buffer to store intermediate coordinates
193  * \param[in]  upd               Update object.
194  * \param[in]  constr            Constraints object.
195  * \param[in]  nullSignaller     Simulation signaller.
196  * \param[in]  trotter_seq       NPT variables.
197  * \param[in]  nrnb              Cycle counters.
198  * \param[in]  wcycle            Wall-clock cycle counter.
199  */
200 void integrateVVSecondStep(int64_t                                  step,
201                            const t_inputrec*                        ir,
202                            t_forcerec*                              fr,
203                            t_commrec*                               cr,
204                            t_state*                                 state,
205                            t_mdatoms*                               mdatoms,
206                            const t_fcdata&                          fcdata,
207                            t_extmass*                               MassQ,
208                            t_vcm*                                   vcm,
209                            pull_t*                                  pull_work,
210                            gmx_enerdata_t*                          enerd,
211                            gmx_ekindata_t*                          ekind,
212                            gmx_global_stat*                         gstat,
213                            real*                                    dvdl_constr,
214                            bool                                     bCalcVir,
215                            tensor                                   total_vir,
216                            tensor                                   shake_vir,
217                            tensor                                   force_vir,
218                            tensor                                   pres,
219                            matrix                                   M,
220                            matrix                                   lastbox,
221                            bool                                     do_log,
222                            bool                                     do_ene,
223                            bool                                     bGStat,
224                            bool*                                    bSumEkinhOld,
225                            gmx::ForceBuffers*                       f,
226                            std::vector<gmx::RVec>*                  cbuf,
227                            gmx::Update*                             upd,
228                            gmx::Constraints*                        constr,
229                            gmx::SimulationSignaller*                nullSignaller,
230                            std::array<std::vector<int>, ettTSEQMAX> trotter_seq,
231                            t_nrnb*                                  nrnb,
232                            gmx_wallcycle*                           wcycle);
233
234
235 #endif // GMX_MDLIB_UPDATE_VV_H