2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
38 #ifndef GMX_MDLIB_UPDATE_H
39 #define GMX_MDLIB_UPDATE_H
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"
52 struct gmx_enerdata_t;
60 enum class ParticleType;
69 * \brief Contains data for update phase */
73 /*! \brief Constructor
75 * \param[in] inputRecord Input record, used to construct SD object.
76 * \param[in] boxDeformation Periodic box deformation object.
78 Update(const t_inputrec& inputRecord, BoxDeformation* boxDeformation);
81 /*! \brief Get the pointer to updated coordinates
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.
87 * \returns The pointer to the intermediate coordinates buffer.
89 PaddedVector<gmx::RVec>* xp();
90 /*!\brief Getter to local copy of box deformation class.
92 * \returns handle to box deformation class
94 BoxDeformation* deform() const;
95 /*! \brief Sets data that changes only at domain decomposition time.
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
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);
107 /*! \brief Perform numerical integration step.
109 * Selects the appropriate integrator, based on the input record and performs a numerical integration step.
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.
127 void update_coords(const t_inputrec& inputRecord,
130 bool havePartiallyFrozenAtoms,
131 gmx::ArrayRef<const ParticleType> ptype,
132 gmx::ArrayRef<const real> invMass,
133 gmx::ArrayRef<const rvec> invMassPerDim,
135 const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
137 const gmx_ekindata_t* ekind,
141 bool haveConstraints);
143 /*! \brief Finalize the coordinate update.
145 * Copy the updated coordinates to the main coordinates buffer for the atoms that are not frozen.
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.
154 void finish_update(const t_inputrec& inputRecord,
155 bool havePartiallyFrozenAtoms,
158 gmx_wallcycle* wcycle,
159 bool haveConstraints);
161 /*! \brief Secong part of the SD integrator.
163 * The first part of integration is performed in the update_coords(...) method.
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.
181 void update_sd_second_half(const t_inputrec& inputRecord,
185 gmx::ArrayRef<const ParticleType> ptype,
186 gmx::ArrayRef<const real> invMass,
190 gmx_wallcycle* wcycle,
191 gmx::Constraints* constr,
195 /*! \brief Performs a leap-frog update without updating \p state so the constrain virial
198 void update_for_constraint_virial(const t_inputrec& inputRecord,
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);
207 /*! \brief Update pre-computed constants that depend on the reference temperature for coupling.
209 * This could change e.g. in simulated annealing.
211 * \param[in] inputRecord Input record.
213 void update_temperature_constants(const t_inputrec& inputRecord);
215 /*!\brief Getter for the list of the randomize groups.
217 * Needed for Andersen temperature control.
219 * \returns Reference to the groups from the SD data object.
221 const std::vector<bool>& getAndersenRandomizeGroup() const;
222 /*!\brief Getter for the list of the Boltzmann factors.
224 * Needed for Andersen temperature control.
226 * \returns Reference to the Boltzmann factors from the SD data object.
228 const std::vector<real>& getBoltzmanFactor() const;
231 //! Implementation type.
233 //! Implementation object.
234 std::unique_ptr<Impl> impl_;
240 * Compute the partial kinetic energy for home particles;
241 * will be accumulated in the calling routine.
244 * Ekin = SUM(i) 0.5 m[i] v[i] (x) v[i]
246 * use v[i] = v[i] - u[i] when calculating temperature
248 * u must be accumulated already.
250 * Now also computes the contribution of the kinetic energy to the
256 void init_ekinstate(ekinstate_t* ekinstate, const t_inputrec* ir);
258 void update_ekinstate(ekinstate_t* ekinstate, const gmx_ekindata_t* ekind);
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);
264 /*! \brief Computes the atom range for a thread to operate on, ensuring SIMD aligned ranges
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
272 void getThreadAtomRange(int numThreads, int threadIndex, int numAtoms, int* startAtom, int* endAtom);