3f122b59f7b1ebac6f63e5e5e6266f6857885325
[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, 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 "gromacs/math/paddedvector.h"
42 #include "gromacs/math/vectypes.h"
43 #include "gromacs/mdtypes/md_enums.h"
44 #include "gromacs/timing/wallcycle.h"
45 #include "gromacs/utility/arrayref.h"
46 #include "gromacs/utility/basedefinitions.h"
47 #include "gromacs/utility/classhelpers.h"
48 #include "gromacs/utility/real.h"
49
50 class ekinstate_t;
51 struct 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_mdatoms;
59 struct t_nrnb;
60 class t_state;
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 Resizes buffer that stores intermediate coordinates.
96      *
97      * \param[in] numAtoms  Updated number of atoms.
98      */
99     void setNumAtoms(int numAtoms);
100
101     /*! \brief Perform numerical integration step.
102      *
103      * Selects the appropriate integrator, based on the input record and performs a numerical integration step.
104      *
105      * \param[in]  inputRecord      Input record.
106      * \param[in]  step             Current timestep.
107      * \param[in]  md               MD atoms data.
108      * \param[in]  state            System state object.
109      * \param[in]  f                Buffer with atomic forces for home particles.
110      * \param[in]  fcdata           Force calculation data to update distance and orientation restraints.
111      * \param[in]  ekind            Kinetic energy data (for temperature coupling, energy groups, etc.).
112      * \param[in]  M                Parrinello-Rahman velocity scaling matrix.
113      * \param[in]  updatePart       What should be updated, coordinates or velocities. This enum only used in VV integrator.
114      * \param[in]  cr               Comunication record  (Old comment: these shouldn't be here -- need to think about it).
115      * \param[in]  haveConstraints  If the system has constraints.
116      */
117     void update_coords(const t_inputrec&                                inputRecord,
118                        int64_t                                          step,
119                        const t_mdatoms*                                 md,
120                        t_state*                                         state,
121                        const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
122                        const t_fcdata&                                  fcdata,
123                        const gmx_ekindata_t*                            ekind,
124                        const matrix                                     M,
125                        int                                              updatePart,
126                        const t_commrec*                                 cr,
127                        bool                                             haveConstraints);
128
129     /*! \brief Finalize the coordinate update.
130      *
131      * Copy the updated coordinates to the main coordinates buffer for the atoms that are not frozen.
132      *
133      * \param[in]  inputRecord      Input record.
134      * \param[in]  md               MD atoms data.
135      * \param[in]  state            System state object.
136      * \param[in]  wcycle           Wall-clock cycle counter.
137      * \param[in]  haveConstraints  If the system has constraints.
138      */
139     void finish_update(const t_inputrec& inputRecord,
140                        const t_mdatoms*  md,
141                        t_state*          state,
142                        gmx_wallcycle_t   wcycle,
143                        bool              haveConstraints);
144
145     /*! \brief Secong part of the SD integrator.
146      *
147      * The first part of integration is performed in the update_coords(...) method.
148      *
149      * \param[in]  inputRecord  Input record.
150      * \param[in]  step         Current timestep.
151      * \param[in]  dvdlambda    Free energy derivative. Contribution to be added to the bonded
152      * interactions. \param[in]  md           MD atoms data. \param[in]  state        System state
153      * object. \param[in]  cr           Comunication record. \param[in]  nrnb         Cycle
154      * counters. \param[in]  wcycle       Wall-clock cycle counter. \param[in]  constr Constraints
155      * object. The constraints are applied on coordinates after update. \param[in]  do_log       If
156      * this is logging step. \param[in]  do_ene       If this is an energy evaluation step.
157      */
158     void update_sd_second_half(const t_inputrec& inputRecord,
159                                int64_t           step,
160                                real*             dvdlambda,
161                                const t_mdatoms*  md,
162                                t_state*          state,
163                                const t_commrec*  cr,
164                                t_nrnb*           nrnb,
165                                gmx_wallcycle_t   wcycle,
166                                gmx::Constraints* constr,
167                                bool              do_log,
168                                bool              do_ene);
169     /*! \brief Update pre-computed constants that depend on the reference temperature for coupling.
170      *
171      * This could change e.g. in simulated annealing.
172      *
173      * \param[in]  inputRecord  Input record.
174      */
175     void update_temperature_constants(const t_inputrec& inputRecord);
176
177     /*!\brief Getter for the list of the randomize groups.
178      *
179      *  Needed for Andersen temperature control.
180      *
181      * \returns Reference to the groups from the SD data object.
182      */
183     const std::vector<bool>& getAndersenRandomizeGroup() const;
184     /*!\brief Getter for the list of the Boltzmann factors.
185      *
186      *  Needed for Andersen temperature control.
187      *
188      * \returns Reference to the Boltzmann factors from the SD data object.
189      */
190     const std::vector<real>& getBoltzmanFactor() const;
191
192 private:
193     //! Implementation type.
194     class Impl;
195     //! Implementation object.
196     PrivateImplPointer<Impl> impl_;
197 };
198
199 }; // namespace gmx
200
201 /*
202  * Compute the partial kinetic energy for home particles;
203  * will be accumulated in the calling routine.
204  * The tensor is
205  *
206  * Ekin = SUM(i) 0.5 m[i] v[i] (x) v[i]
207  *
208  *     use v[i] = v[i] - u[i] when calculating temperature
209  *
210  * u must be accumulated already.
211  *
212  * Now also computes the contribution of the kinetic energy to the
213  * free energy
214  *
215  */
216
217
218 void init_ekinstate(ekinstate_t* ekinstate, const t_inputrec* ir);
219
220 void update_ekinstate(ekinstate_t* ekinstate, const gmx_ekindata_t* ekind);
221
222 /*! \brief Restores data from \p ekinstate to \p ekind, then broadcasts it
223    to the rest of the simulation */
224 void restore_ekinstate_from_state(const t_commrec* cr, gmx_ekindata_t* ekind, const ekinstate_t* ekinstate);
225
226 /*! \brief Computes the atom range for a thread to operate on, ensuring SIMD aligned ranges
227  *
228  * \param[in]  numThreads   The number of threads to divide atoms over
229  * \param[in]  threadIndex  The thread to get the range for
230  * \param[in]  numAtoms     The total number of atoms (on this rank)
231  * \param[out] startAtom    The start of the atom range
232  * \param[out] endAtom      The end of the atom range, note that this is in general not a multiple of the SIMD width
233  */
234 void getThreadAtomRange(int numThreads, int threadIndex, int numAtoms, int* startAtom, int* endAtom);
235
236 #endif