Reimplement constant acceleration groups
[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 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_H
39 #define GMX_MDLIB_UPDATE_H
40
41 #include <memory>
42
43 #include "gromacs/math/paddedvector.h"
44 #include "gromacs/math/vectypes.h"
45 #include "gromacs/mdtypes/md_enums.h"
46 #include "gromacs/timing/wallcycle.h"
47 #include "gromacs/utility/arrayref.h"
48 #include "gromacs/utility/real.h"
49
50 class ekinstate_t;
51 class gmx_ekindata_t;
52 struct gmx_enerdata_t;
53 enum class PbcType;
54 struct t_fcdata;
55 struct t_graph;
56 struct t_grpopts;
57 struct t_inputrec;
58 struct t_nrnb;
59 class t_state;
60 enum class ParticleType;
61
62 namespace gmx
63 {
64 class BoxDeformation;
65 class Constraints;
66
67
68 /*! \libinternal
69  * \brief Contains data for update phase */
70 class Update
71 {
72 public:
73     /*! \brief Constructor
74      *
75      * \param[in] inputRecord     Input record, used to construct SD object.
76      * \param[in] boxDeformation  Periodic box deformation object.
77      */
78     Update(const t_inputrec& inputRecord, BoxDeformation* boxDeformation);
79     //! Destructor
80     ~Update();
81     /*! \brief Get the pointer to updated coordinates
82      *
83      * Update saves the updated coordinates into separate buffer, so that constraints will have
84      * access to both updated and not update coordinates. For that, update owns a separate buffer.
85      * See finish_update(...) for details.
86      *
87      * \returns The pointer to the intermediate coordinates buffer.
88      */
89     PaddedVector<gmx::RVec>* xp();
90     /*!\brief Getter to local copy of box deformation class.
91      *
92      * \returns handle to box deformation class
93      */
94     BoxDeformation* deform() const;
95     /*! \brief Sets data that changes only at domain decomposition time.
96      *
97      * \param[in] numAtoms  Updated number of atoms.
98      * \param[in] cFREEZE   Group index for freezing
99      * \param[in] cTC       Group index for center of mass motion removal
100      * \param[in] cAcceleration  Group index for constant acceleration groups
101      */
102     void updateAfterPartition(int                                 numAtoms,
103                               gmx::ArrayRef<const unsigned short> cFREEZE,
104                               gmx::ArrayRef<const unsigned short> cTC,
105                               gmx::ArrayRef<const unsigned short> cAcceleration);
106
107     /*! \brief Perform numerical integration step.
108      *
109      * Selects the appropriate integrator, based on the input record and performs a numerical integration step.
110      *
111      * \param[in]  inputRecord               Input record.
112      * \param[in]  step                      Current timestep.
113      * \param[in]  homenr                    The number of atoms on this processor.
114      * \param[in]  havePartiallyFrozenAtoms  Whether atoms are frozen along 1 or 2 (not 3) dimensions?
115      * \param[in]  ptype                     The list of particle types.
116      * \param[in]  invMass                   Inverse atomic mass per atom, 0 for vsites and shells.
117      * \param[in]  invMassPerDim             Inverse atomic mass per atom and dimension, 0 for vsites, shells and frozen dimensions
118      * \param[in]  state                     System state object.
119      * \param[in]  f                         Buffer with atomic forces for home particles.
120      * \param[in]  fcdata                    Force calculation data to update distance and orientation restraints.
121      * \param[in]  ekind                     Kinetic energy data (for temperature coupling, energy groups, etc.).
122      * \param[in]  M                         Parrinello-Rahman velocity scaling matrix.
123      * \param[in]  updatePart                What should be updated, coordinates or velocities. This enum only used in VV integrator.
124      * \param[in]  cr                        Comunication record  (Old comment: these shouldn't be here -- need to think about it).
125      * \param[in]  haveConstraints           If the system has constraints.
126      */
127     void update_coords(const t_inputrec&                                inputRecord,
128                        int64_t                                          step,
129                        int                                              homenr,
130                        bool                                             havePartiallyFrozenAtoms,
131                        gmx::ArrayRef<const ParticleType>                ptype,
132                        gmx::ArrayRef<const real>                        invMass,
133                        gmx::ArrayRef<const rvec>                        invMassPerDim,
134                        t_state*                                         state,
135                        const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
136                        t_fcdata*                                        fcdata,
137                        const gmx_ekindata_t*                            ekind,
138                        const matrix                                     M,
139                        int                                              updatePart,
140                        const t_commrec*                                 cr,
141                        bool                                             haveConstraints);
142
143     /*! \brief Finalize the coordinate update.
144      *
145      * Copy the updated coordinates to the main coordinates buffer for the atoms that are not frozen.
146      *
147      * \param[in]  inputRecord      Input record.
148      * \param[in]  havePartiallyFrozenAtoms  Whether atoms are frozen along 1 or 2 (not 3) dimensions?
149      * \param[in]  homenr                    The number of atoms on this processor.
150      * \param[in]  state            System state object.
151      * \param[in]  wcycle           Wall-clock cycle counter.
152      * \param[in]  haveConstraints  If the system has constraints.
153      */
154     void finish_update(const t_inputrec& inputRecord,
155                        bool              havePartiallyFrozenAtoms,
156                        int               homenr,
157                        t_state*          state,
158                        gmx_wallcycle*    wcycle,
159                        bool              haveConstraints);
160
161     /*! \brief Secong part of the SD integrator.
162      *
163      * The first part of integration is performed in the update_coords(...) method.
164      *
165      * \param[in]  inputRecord  Input record.
166      * \param[in]  step         Current timestep.
167      * \param[in]  dvdlambda    Free energy derivative. Contribution to be added to
168      *                          the bonded interactions.
169      * \param[in]  homenr       The number of atoms on this processor.
170      * \param[in]  ptype        The list of particle types.
171      * \param[in]  invMass      Inverse atomic mass per atom, 0 for vsites and shells.
172      * \param[in]  state        System state object.
173      * \param[in]  cr           Comunication record.
174      * \param[in]  nrnb         Cycle counters.
175      * \param[in]  wcycle       Wall-clock cycle counter.
176      * \param[in]  constr       Constraints object. The constraints are applied
177      *                          on coordinates after update.
178      * \param[in]  do_log       If this is logging step.
179      * \param[in]  do_ene       If this is an energy evaluation step.
180      */
181     void update_sd_second_half(const t_inputrec&                 inputRecord,
182                                int64_t                           step,
183                                real*                             dvdlambda,
184                                int                               homenr,
185                                gmx::ArrayRef<const ParticleType> ptype,
186                                gmx::ArrayRef<const real>         invMass,
187                                t_state*                          state,
188                                const t_commrec*                  cr,
189                                t_nrnb*                           nrnb,
190                                gmx_wallcycle*                    wcycle,
191                                gmx::Constraints*                 constr,
192                                bool                              do_log,
193                                bool                              do_ene);
194
195     /*! \brief Performs a leap-frog update without updating \p state so the constrain virial
196      * can be computed.
197      */
198     void update_for_constraint_virial(const t_inputrec&         inputRecord,
199                                       int                       homenr,
200                                       bool                      havePartiallyFrozenAtoms,
201                                       gmx::ArrayRef<const real> invmass,
202                                       gmx::ArrayRef<const rvec> invMassPerDim,
203                                       const t_state&            state,
204                                       const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
205                                       const gmx_ekindata_t&                            ekind);
206
207     /*! \brief Update pre-computed constants that depend on the reference temperature for coupling.
208      *
209      * This could change e.g. in simulated annealing.
210      *
211      * \param[in]  inputRecord  Input record.
212      */
213     void update_temperature_constants(const t_inputrec& inputRecord);
214
215     /*!\brief Getter for the list of the randomize groups.
216      *
217      *  Needed for Andersen temperature control.
218      *
219      * \returns Reference to the groups from the SD data object.
220      */
221     const std::vector<bool>& getAndersenRandomizeGroup() const;
222     /*!\brief Getter for the list of the Boltzmann factors.
223      *
224      *  Needed for Andersen temperature control.
225      *
226      * \returns Reference to the Boltzmann factors from the SD data object.
227      */
228     const std::vector<real>& getBoltzmanFactor() const;
229
230 private:
231     //! Implementation type.
232     class Impl;
233     //! Implementation object.
234     std::unique_ptr<Impl> impl_;
235 };
236
237 }; // namespace gmx
238
239 /*
240  * Compute the partial kinetic energy for home particles;
241  * will be accumulated in the calling routine.
242  * The tensor is
243  *
244  * Ekin = SUM(i) 0.5 m[i] v[i] (x) v[i]
245  *
246  *     use v[i] = v[i] - u[i] when calculating temperature
247  *
248  * u must be accumulated already.
249  *
250  * Now also computes the contribution of the kinetic energy to the
251  * free energy
252  *
253  */
254
255
256 void init_ekinstate(ekinstate_t* ekinstate, const t_inputrec* ir);
257
258 void update_ekinstate(ekinstate_t* ekinstate, const gmx_ekindata_t* ekind);
259
260 /*! \brief Restores data from \p ekinstate to \p ekind, then broadcasts it
261    to the rest of the simulation */
262 void restore_ekinstate_from_state(const t_commrec* cr, gmx_ekindata_t* ekind, const ekinstate_t* ekinstate);
263
264 /*! \brief Computes the atom range for a thread to operate on, ensuring SIMD aligned ranges
265  *
266  * \param[in]  numThreads   The number of threads to divide atoms over
267  * \param[in]  threadIndex  The thread to get the range for
268  * \param[in]  numAtoms     The total number of atoms (on this rank)
269  * \param[out] startAtom    The start of the atom range
270  * \param[out] endAtom      The end of the atom range, note that this is in general not a multiple of the SIMD width
271  */
272 void getThreadAtomRange(int numThreads, int threadIndex, int numAtoms, int* startAtom, int* endAtom);
273
274 #endif