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