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.
48 #include "gromacs/domdec/domdec_struct.h"
49 #include "gromacs/fileio/confio.h"
50 #include "gromacs/gmxlib/network.h"
51 #include "gromacs/gmxlib/nrnb.h"
52 #include "gromacs/listed_forces/disre.h"
53 #include "gromacs/listed_forces/orires.h"
54 #include "gromacs/math/functions.h"
55 #include "gromacs/math/invertmatrix.h"
56 #include "gromacs/math/paddedvector.h"
57 #include "gromacs/math/units.h"
58 #include "gromacs/math/vec.h"
59 #include "gromacs/math/vecdump.h"
60 #include "gromacs/mdlib/boxdeformation.h"
61 #include "gromacs/mdlib/constr.h"
62 #include "gromacs/mdlib/gmx_omp_nthreads.h"
63 #include "gromacs/mdlib/mdatoms.h"
64 #include "gromacs/mdlib/stat.h"
65 #include "gromacs/mdlib/tgroup.h"
66 #include "gromacs/mdtypes/commrec.h"
67 #include "gromacs/mdtypes/fcdata.h"
68 #include "gromacs/mdtypes/group.h"
69 #include "gromacs/mdtypes/inputrec.h"
70 #include "gromacs/mdtypes/md_enums.h"
71 #include "gromacs/mdtypes/mdatom.h"
72 #include "gromacs/mdtypes/state.h"
73 #include "gromacs/pbcutil/boxutilities.h"
74 #include "gromacs/pbcutil/pbc.h"
75 #include "gromacs/pulling/pull.h"
76 #include "gromacs/random/tabulatednormaldistribution.h"
77 #include "gromacs/random/threefry.h"
78 #include "gromacs/simd/simd.h"
79 #include "gromacs/timing/wallcycle.h"
80 #include "gromacs/topology/atoms.h"
81 #include "gromacs/utility/exceptions.h"
82 #include "gromacs/utility/fatalerror.h"
83 #include "gromacs/utility/futil.h"
84 #include "gromacs/utility/gmxassert.h"
85 #include "gromacs/utility/gmxomp.h"
86 #include "gromacs/utility/smalloc.h"
88 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
103 std::vector<real> bd_rf;
105 std::vector<gmx_sd_const_t> sdc;
106 std::vector<gmx_sd_sigma_t> sdsig;
107 /* andersen temperature control stuff */
108 std::vector<bool> randomize_group;
109 std::vector<real> boltzfac;
111 explicit gmx_stochd_t(const t_inputrec& inputRecord);
114 //! pImpled implementation for Update
119 Impl(const t_inputrec& inputRecord, BoxDeformation* boxDeformation);
123 void update_coords(const t_inputrec& inputRecord,
127 const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
128 const t_fcdata& fcdata,
129 const gmx_ekindata_t* ekind,
133 bool haveConstraints);
135 void finish_update(const t_inputrec& inputRecord,
138 gmx_wallcycle_t wcycle,
139 bool haveConstraints);
141 void update_sd_second_half(const t_inputrec& inputRecord,
148 gmx_wallcycle_t wcycle,
149 gmx::Constraints* constr,
153 void update_for_constraint_virial(const t_inputrec& inputRecord,
155 const t_state& state,
156 const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
157 const gmx_ekindata_t& ekind);
159 void update_temperature_constants(const t_inputrec& inputRecord);
161 const std::vector<bool>& getAndersenRandomizeGroup() const { return sd_.randomize_group; }
163 const std::vector<real>& getBoltzmanFactor() const { return sd_.boltzfac; }
165 PaddedVector<RVec>* xp() { return &xp_; }
167 BoxDeformation* deform() const { return deform_; }
170 //! stochastic dynamics struct
172 //! xprime for constraint algorithms
173 PaddedVector<RVec> xp_;
174 //! Box deformation handler (or nullptr if inactive).
175 BoxDeformation* deform_ = nullptr;
178 Update::Update(const t_inputrec& inputRecord, BoxDeformation* boxDeformation) :
179 impl_(new Impl(inputRecord, boxDeformation)){};
183 const std::vector<bool>& Update::getAndersenRandomizeGroup() const
185 return impl_->getAndersenRandomizeGroup();
188 const std::vector<real>& Update::getBoltzmanFactor() const
190 return impl_->getBoltzmanFactor();
193 PaddedVector<RVec>* Update::xp()
198 BoxDeformation* Update::deform() const
200 return impl_->deform();
203 void Update::update_coords(const t_inputrec& inputRecord,
207 const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
208 const t_fcdata& fcdata,
209 const gmx_ekindata_t* ekind,
213 const bool haveConstraints)
215 return impl_->update_coords(
216 inputRecord, step, md, state, f, fcdata, ekind, M, updatePart, cr, haveConstraints);
219 void Update::finish_update(const t_inputrec& inputRecord,
222 gmx_wallcycle_t wcycle,
223 const bool haveConstraints)
225 return impl_->finish_update(inputRecord, md, state, wcycle, haveConstraints);
228 void Update::update_sd_second_half(const t_inputrec& inputRecord,
235 gmx_wallcycle_t wcycle,
236 gmx::Constraints* constr,
240 return impl_->update_sd_second_half(
241 inputRecord, step, dvdlambda, md, state, cr, nrnb, wcycle, constr, do_log, do_ene);
244 void Update::update_for_constraint_virial(const t_inputrec& inputRecord,
246 const t_state& state,
247 const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
248 const gmx_ekindata_t& ekind)
250 return impl_->update_for_constraint_virial(inputRecord, md, state, f, ekind);
253 void Update::update_temperature_constants(const t_inputrec& inputRecord)
255 return impl_->update_temperature_constants(inputRecord);
258 /*! \brief Sets whether we store the updated velocities */
259 enum class StoreUpdatedVelocities
261 yes, //!< Store the updated velocities
262 no //!< Do not store the updated velocities
265 /*! \brief Sets the number of different temperature coupling values */
266 enum class NumTempScaleValues
268 single, //!< Single T-scaling value (either one group or all values =1)
269 multiple //!< Multiple T-scaling values, need to use T-group indices
272 /*! \brief Sets if to apply no or diagonal Parrinello-Rahman pressure scaling
274 * Note that this enum is only used in updateMDLeapfrogSimple(), which does
275 * not handle fully anistropic Parrinello-Rahman scaling, so we only have
276 * options \p no and \p diagonal here and no anistropic option.
278 enum class ApplyParrinelloRahmanVScaling
280 no, //!< Do not apply velocity scaling (not a PR-coupling run or step)
281 diagonal //!< Apply velocity scaling using a diagonal matrix
284 /*! \brief Integrate using leap-frog with T-scaling and optionally diagonal Parrinello-Rahman p-coupling
286 * \tparam storeUpdatedVelocities Tells whether we should store the updated velocities
287 * \tparam numTempScaleValues The number of different T-couple values
288 * \tparam applyPRVScaling Apply Parrinello-Rahman velocity scaling
289 * \param[in] start Index of first atom to update
290 * \param[in] nrend Last atom to update: \p nrend - 1
291 * \param[in] dt The time step
292 * \param[in] dtPressureCouple Time step for pressure coupling
293 * \param[in] invMassPerDim 1/mass per atom and dimension
294 * \param[in] tcstat Temperature coupling information
295 * \param[in] cTC T-coupling group index per atom
296 * \param[in] pRVScaleMatrixDiagonal Parrinello-Rahman v-scale matrix diagonal
297 * \param[in] x Input coordinates
298 * \param[out] xprime Updated coordinates
299 * \param[inout] v Velocities, type either rvec* or const rvec*
300 * \param[in] f Forces
302 * We expect this template to get good SIMD acceleration by most compilers,
303 * unlike the more complex general template.
304 * Note that we might get even better SIMD acceleration when we introduce
305 * aligned (and padded) memory, possibly with some hints for the compilers.
307 template<StoreUpdatedVelocities storeUpdatedVelocities, NumTempScaleValues numTempScaleValues, ApplyParrinelloRahmanVScaling applyPRVScaling, typename VelocityType>
308 static std::enable_if_t<std::is_same<VelocityType, rvec*>::value || std::is_same<VelocityType, const rvec*>::value, void>
309 updateMDLeapfrogSimple(int start,
312 real dtPressureCouple,
313 const rvec* gmx_restrict invMassPerDim,
314 gmx::ArrayRef<const t_grp_tcstat> tcstat,
315 const unsigned short* cTC,
316 const rvec pRVScaleMatrixDiagonal,
317 const rvec* gmx_restrict x,
318 rvec* gmx_restrict xprime,
319 VelocityType gmx_restrict v,
320 const rvec* gmx_restrict f)
324 if (numTempScaleValues == NumTempScaleValues::single)
326 lambdaGroup = tcstat[0].lambda;
329 for (int a = start; a < nrend; a++)
331 if (numTempScaleValues == NumTempScaleValues::multiple)
333 lambdaGroup = tcstat[cTC[a]].lambda;
336 for (int d = 0; d < DIM; d++)
338 /* Note that using rvec invMassPerDim results in more efficient
339 * SIMD code, but this increases the cache pressure.
340 * For large systems with PME on the CPU this slows down the
341 * (then already slow) update by 20%. If all data remains in cache,
342 * using rvec is much faster.
344 real vNew = lambdaGroup * v[a][d] + f[a][d] * invMassPerDim[a][d] * dt;
346 if (applyPRVScaling == ApplyParrinelloRahmanVScaling::diagonal)
348 vNew -= dtPressureCouple * pRVScaleMatrixDiagonal[d] * v[a][d];
350 // TODO: Remove NOLINTs once clang-tidy is updated to v11, it should be able to handle constexpr.
351 if constexpr (storeUpdatedVelocities == StoreUpdatedVelocities::yes) // NOLINT // NOLINTNEXTLINE
355 xprime[a][d] = x[a][d] + vNew * dt; // NOLINT(readability-misleading-indentation)
360 #if GMX_SIMD && GMX_SIMD_HAVE_REAL
361 # define GMX_HAVE_SIMD_UPDATE 1
363 # define GMX_HAVE_SIMD_UPDATE 0
366 #if GMX_HAVE_SIMD_UPDATE
368 /*! \brief Load (aligned) the contents of GMX_SIMD_REAL_WIDTH rvec elements sequentially into 3 SIMD registers
370 * The loaded output is:
371 * \p r0: { r[index][XX], r[index][YY], ... }
373 * \p r2: { ..., r[index+GMX_SIMD_REAL_WIDTH-1][YY], r[index+GMX_SIMD_REAL_WIDTH-1][ZZ] }
375 * \param[in] r Real to an rvec array, has to be aligned to SIMD register width
376 * \param[in] index Index of the first rvec triplet of reals to load
377 * \param[out] r0 Pointer to first SIMD register
378 * \param[out] r1 Pointer to second SIMD register
379 * \param[out] r2 Pointer to third SIMD register
381 static inline void simdLoadRvecs(const rvec* r, int index, SimdReal* r0, SimdReal* r1, SimdReal* r2)
383 const real* realPtr = r[index];
385 GMX_ASSERT(isSimdAligned(realPtr), "Pointer should be SIMD aligned");
387 *r0 = simdLoad(realPtr + 0 * GMX_SIMD_REAL_WIDTH);
388 *r1 = simdLoad(realPtr + 1 * GMX_SIMD_REAL_WIDTH);
389 *r2 = simdLoad(realPtr + 2 * GMX_SIMD_REAL_WIDTH);
392 /*! \brief Store (aligned) 3 SIMD registers sequentially to GMX_SIMD_REAL_WIDTH rvec elements
394 * The stored output is:
395 * \p r[index] = { { r0[0], r0[1], ... }
397 * \p r[index+GMX_SIMD_REAL_WIDTH-1] = { ... , r2[GMX_SIMD_REAL_WIDTH-2], r2[GMX_SIMD_REAL_WIDTH-1] }
399 * \param[out] r Pointer to an rvec array, has to be aligned to SIMD register width
400 * \param[in] index Index of the first rvec triplet of reals to store to
401 * \param[in] r0 First SIMD register
402 * \param[in] r1 Second SIMD register
403 * \param[in] r2 Third SIMD register
405 static inline void simdStoreRvecs(rvec* r, int index, SimdReal r0, SimdReal r1, SimdReal r2)
407 real* realPtr = r[index];
409 GMX_ASSERT(isSimdAligned(realPtr), "Pointer should be SIMD aligned");
411 store(realPtr + 0 * GMX_SIMD_REAL_WIDTH, r0);
412 store(realPtr + 1 * GMX_SIMD_REAL_WIDTH, r1);
413 store(realPtr + 2 * GMX_SIMD_REAL_WIDTH, r2);
416 /*! \brief Integrate using leap-frog with single group T-scaling and SIMD
418 * \tparam storeUpdatedVelocities Tells whether we should store the updated velocities
419 * \param[in] start Index of first atom to update
420 * \param[in] nrend Last atom to update: \p nrend - 1
421 * \param[in] dt The time step
422 * \param[in] invMass 1/mass per atom
423 * \param[in] tcstat Temperature coupling information
424 * \param[in] x Input coordinates
425 * \param[out] xprime Updated coordinates
426 * \param[inout] v Velocities, type either rvec* or const rvec*
427 * \param[in] f Forces
429 template<StoreUpdatedVelocities storeUpdatedVelocities, typename VelocityType>
430 static std::enable_if_t<std::is_same<VelocityType, rvec*>::value || std::is_same<VelocityType, const rvec*>::value, void>
431 updateMDLeapfrogSimpleSimd(int start,
434 const real* gmx_restrict invMass,
435 gmx::ArrayRef<const t_grp_tcstat> tcstat,
436 const rvec* gmx_restrict x,
437 rvec* gmx_restrict xprime,
438 VelocityType gmx_restrict v,
439 const rvec* gmx_restrict f)
441 SimdReal timestep(dt);
442 SimdReal lambdaSystem(tcstat[0].lambda);
444 /* We declare variables here, since code is often slower when declaring them inside the loop */
446 /* Note: We should implement a proper PaddedVector, so we don't need this check */
447 GMX_ASSERT(isSimdAligned(invMass), "invMass should be aligned");
449 for (int a = start; a < nrend; a += GMX_SIMD_REAL_WIDTH)
451 SimdReal invMass0, invMass1, invMass2;
452 expandScalarsToTriplets(simdLoad(invMass + a), &invMass0, &invMass1, &invMass2);
456 simdLoadRvecs(v, a, &v0, &v1, &v2);
457 simdLoadRvecs(f, a, &f0, &f1, &f2);
459 v0 = fma(f0 * invMass0, timestep, lambdaSystem * v0);
460 v1 = fma(f1 * invMass1, timestep, lambdaSystem * v1);
461 v2 = fma(f2 * invMass2, timestep, lambdaSystem * v2);
463 // TODO: Remove NOLINTs once clang-tidy is updated to v11, it should be able to handle constexpr.
464 if constexpr (storeUpdatedVelocities == StoreUpdatedVelocities::yes) // NOLINT // NOLINTNEXTLINE
466 simdStoreRvecs(v, a, v0, v1, v2);
469 SimdReal x0, x1, x2; // NOLINT(readability-misleading-indentation)
470 simdLoadRvecs(x, a, &x0, &x1, &x2);
472 SimdReal xprime0 = fma(v0, timestep, x0);
473 SimdReal xprime1 = fma(v1, timestep, x1);
474 SimdReal xprime2 = fma(v2, timestep, x2);
476 simdStoreRvecs(xprime, a, xprime0, xprime1, xprime2);
480 #endif // GMX_HAVE_SIMD_UPDATE
482 /*! \brief Sets the NEMD acceleration type */
483 enum class AccelerationType
489 /*! \brief Integrate using leap-frog with support for everything.
491 * \tparam accelerationType Type of NEMD acceleration.
492 * \param[in] start Index of first atom to update.
493 * \param[in] nrend Last atom to update: \p nrend - 1.
494 * \param[in] doNoseHoover If to apply Nose-Hoover T-coupling.
495 * \param[in] dt The time step.
496 * \param[in] dtPressureCouple Time step for pressure coupling, is 0 when no pressure
497 * coupling should be applied at this step.
498 * \param[in] md Atom properties.
499 * \param[in] ekind Kinetic energy data.
500 * \param[in] box The box dimensions.
501 * \param[in] x Input coordinates.
502 * \param[out] xprime Updated coordinates.
503 * \param[inout] v Velocities.
504 * \param[in] f Forces.
505 * \param[in] nh_vxi Nose-Hoover velocity scaling factors.
506 * \param[in] nsttcouple Frequency of the temperature coupling steps.
507 * \param[in] M Parrinello-Rahman scaling matrix.
509 template<AccelerationType accelerationType>
510 static void updateMDLeapfrogGeneral(int start,
514 real dtPressureCouple,
516 const gmx_ekindata_t* ekind,
518 const rvec* gmx_restrict x,
519 rvec* gmx_restrict xprime,
520 rvec* gmx_restrict v,
521 const rvec* gmx_restrict f,
522 const double* gmx_restrict nh_vxi,
523 const int nsttcouple,
526 /* This is a version of the leap-frog integrator that supports
527 * all combinations of T-coupling, P-coupling and NEMD.
528 * Nose-Hoover uses the reversible leap-frog integrator from
529 * Holian et al. Phys Rev E 52(3) : 2338, 1995
532 gmx::ArrayRef<const t_grp_tcstat> tcstat = ekind->tcstat;
533 const unsigned short* cTC = md->cTC;
535 const rvec* gmx_restrict invMassPerDim = md->invMassPerDim;
537 /* Initialize group values, changed later when multiple groups are used */
540 real omega_Z = 2 * static_cast<real>(M_PI) / box[ZZ][ZZ];
542 for (int n = start; n < nrend; n++)
548 real lg = tcstat[gt].lambda;
551 real cosineZ, vCosine;
552 switch (accelerationType)
554 case AccelerationType::none: copy_rvec(v[n], vRel); break;
555 case AccelerationType::cosine:
556 cosineZ = std::cos(x[n][ZZ] * omega_Z);
557 vCosine = cosineZ * ekind->cosacc.vcos;
558 /* Avoid scaling the cosine profile velocity */
559 copy_rvec(v[n], vRel);
567 /* Here we account for multiple time stepping, by increasing
568 * the Nose-Hoover correction by nsttcouple
569 * TODO: This can be pre-computed.
571 factorNH = 0.5 * nsttcouple * dt * nh_vxi[gt];
574 for (int d = 0; d < DIM; d++)
576 real vNew = (lg * vRel[d]
577 + (f[n][d] * invMassPerDim[n][d] * dt - factorNH * vRel[d]
578 - dtPressureCouple * iprod(M[d], vRel)))
580 switch (accelerationType)
582 case AccelerationType::none: break;
583 case AccelerationType::cosine:
586 /* Add back the mean velocity and apply acceleration */
587 vNew += vCosine + cosineZ * ekind->cosacc.cos_accel * dt;
592 xprime[n][d] = x[n][d] + vNew * dt;
597 /*! \brief Handles the Leap-frog MD x and v integration */
598 static void do_update_md(int start,
602 const rvec* gmx_restrict x,
603 rvec* gmx_restrict xprime,
604 rvec* gmx_restrict v,
605 const rvec* gmx_restrict f,
606 const TemperatureCoupling etc,
607 const PressureCoupling epc,
608 const int nsttcouple,
609 const int nstpcouple,
611 const gmx_ekindata_t* ekind,
613 const double* gmx_restrict nh_vxi,
616 GMX_ASSERT(nrend == start || xprime != x,
617 "For SIMD optimization certain compilers need to have xprime != x");
619 /* Note: Berendsen pressure scaling is handled after do_update_md() */
621 (etc != TemperatureCoupling::No && do_per_step(step + nsttcouple - 1, nsttcouple));
622 bool doNoseHoover = (etc == TemperatureCoupling::NoseHoover && doTempCouple);
623 bool doParrinelloRahman = (epc == PressureCoupling::ParrinelloRahman
624 && do_per_step(step + nstpcouple - 1, nstpcouple));
625 bool doPROffDiagonal = (doParrinelloRahman && (M[YY][XX] != 0 || M[ZZ][XX] != 0 || M[ZZ][YY] != 0));
627 real dtPressureCouple = (doParrinelloRahman ? nstpcouple * dt : 0);
629 /* NEMD cosine acceleration is applied in updateMDLeapFrogGeneral */
630 bool doAcceleration = (ekind->cosacc.cos_accel != 0);
632 if (doNoseHoover || doPROffDiagonal || doAcceleration)
635 if (!doParrinelloRahman)
637 /* We should not apply PR scaling at this step */
647 updateMDLeapfrogGeneral<AccelerationType::none>(
648 start, nrend, doNoseHoover, dt, dtPressureCouple, md, ekind, box, x, xprime, v, f, nh_vxi, nsttcouple, stepM);
652 updateMDLeapfrogGeneral<AccelerationType::cosine>(
653 start, nrend, doNoseHoover, dt, dtPressureCouple, md, ekind, box, x, xprime, v, f, nh_vxi, nsttcouple, stepM);
658 /* We determine if we have a single T-coupling lambda value for all
659 * atoms. That allows for better SIMD acceleration in the template.
660 * If we do not do temperature coupling (in the run or this step),
661 * all scaling values are 1, so we effectively have a single value.
663 bool haveSingleTempScaleValue = (!doTempCouple || ekind->ngtc == 1);
665 /* Extract some pointers needed by all cases */
666 const unsigned short* cTC = md->cTC;
667 gmx::ArrayRef<const t_grp_tcstat> tcstat = ekind->tcstat;
668 const rvec* invMassPerDim = md->invMassPerDim;
670 if (doParrinelloRahman)
672 GMX_ASSERT(!doPROffDiagonal,
673 "updateMDLeapfrogSimple only support diagonal Parrinello-Rahman scaling "
677 for (int d = 0; d < DIM; d++)
682 if (haveSingleTempScaleValue)
684 updateMDLeapfrogSimple<StoreUpdatedVelocities::yes, NumTempScaleValues::single, ApplyParrinelloRahmanVScaling::diagonal>(
685 start, nrend, dt, dtPressureCouple, invMassPerDim, tcstat, cTC, diagM, x, xprime, v, f);
689 updateMDLeapfrogSimple<StoreUpdatedVelocities::yes, NumTempScaleValues::multiple, ApplyParrinelloRahmanVScaling::diagonal>(
690 start, nrend, dt, dtPressureCouple, invMassPerDim, tcstat, cTC, diagM, x, xprime, v, f);
695 if (haveSingleTempScaleValue)
697 /* Note that modern compilers are pretty good at vectorizing
698 * updateMDLeapfrogSimple(). But the SIMD version will still
699 * be faster because invMass lowers the cache pressure
700 * compared to invMassPerDim.
702 #if GMX_HAVE_SIMD_UPDATE
703 /* Check if we can use invmass instead of invMassPerDim */
704 if (!md->havePartiallyFrozenAtoms)
706 updateMDLeapfrogSimpleSimd<StoreUpdatedVelocities::yes>(
707 start, nrend, dt, md->invmass, tcstat, x, xprime, v, f);
712 updateMDLeapfrogSimple<StoreUpdatedVelocities::yes, NumTempScaleValues::single, ApplyParrinelloRahmanVScaling::no>(
713 start, nrend, dt, dtPressureCouple, invMassPerDim, tcstat, cTC, nullptr, x, xprime, v, f);
718 updateMDLeapfrogSimple<StoreUpdatedVelocities::yes, NumTempScaleValues::multiple, ApplyParrinelloRahmanVScaling::no>(
719 start, nrend, dt, dtPressureCouple, invMassPerDim, tcstat, cTC, nullptr, x, xprime, v, f);
724 /*! \brief Handles the Leap-frog MD x and v integration */
725 static void doUpdateMDDoNotUpdateVelocities(int start,
728 const rvec* gmx_restrict x,
729 rvec* gmx_restrict xprime,
730 const rvec* gmx_restrict v,
731 const rvec* gmx_restrict f,
733 const gmx_ekindata_t& ekind)
735 GMX_ASSERT(nrend == start || xprime != x,
736 "For SIMD optimization certain compilers need to have xprime != x");
738 gmx::ArrayRef<const t_grp_tcstat> tcstat = ekind.tcstat;
740 /* Check if we can use invmass instead of invMassPerDim */
741 #if GMX_HAVE_SIMD_UPDATE
742 if (!md.havePartiallyFrozenAtoms)
744 updateMDLeapfrogSimpleSimd<StoreUpdatedVelocities::no>(
745 start, nrend, dt, md.invmass, tcstat, x, xprime, v, f);
750 updateMDLeapfrogSimple<StoreUpdatedVelocities::no, NumTempScaleValues::single, ApplyParrinelloRahmanVScaling::no>(
751 start, nrend, dt, dt, md.invMassPerDim, tcstat, nullptr, nullptr, x, xprime, v, f);
755 static void do_update_vv_vel(int start,
758 const ivec nFreeze[],
759 const real invmass[],
760 const ParticleType ptype[],
761 const unsigned short cFREEZE[],
774 g = 0.25 * dt * veta * alpha;
776 mv2 = gmx::series_sinhx(g);
783 for (n = start; n < nrend; n++)
785 real w_dt = invmass[n] * dt;
791 for (d = 0; d < DIM; d++)
793 if ((ptype[n] != ParticleType::Shell) && !nFreeze[gf][d])
795 v[n][d] = mv1 * (mv1 * v[n][d] + 0.5 * (w_dt * mv2 * f[n][d]));
803 } /* do_update_vv_vel */
805 static void do_update_vv_pos(int start,
808 const ivec nFreeze[],
809 const ParticleType ptype[],
810 const unsigned short cFREEZE[],
821 /* Would it make more sense if Parrinello-Rahman was put here? */
826 mr2 = gmx::series_sinhx(g);
834 for (n = start; n < nrend; n++)
842 for (d = 0; d < DIM; d++)
844 if ((ptype[n] != ParticleType::Shell) && !nFreeze[gf][d])
846 xprime[n][d] = mr1 * (mr1 * x[n][d] + mr2 * dt * v[n][d]);
850 xprime[n][d] = x[n][d];
854 } /* do_update_vv_pos */
856 gmx_stochd_t::gmx_stochd_t(const t_inputrec& inputRecord)
858 const t_grpopts* opts = &inputRecord.opts;
859 const int ngtc = opts->ngtc;
861 if (inputRecord.eI == IntegrationAlgorithm::BD)
865 else if (EI_SD(inputRecord.eI))
870 for (int gt = 0; gt < ngtc; gt++)
872 if (opts->tau_t[gt] > 0)
874 sdc[gt].em = std::exp(-inputRecord.delta_t / opts->tau_t[gt]);
878 /* No friction and noise on this group */
883 else if (ETC_ANDERSEN(inputRecord.etc))
885 randomize_group.resize(ngtc);
886 boltzfac.resize(ngtc);
888 /* for now, assume that all groups, if randomized, are randomized at the same rate, i.e. tau_t is the same. */
889 /* since constraint groups don't necessarily match up with temperature groups! This is checked in readir.c */
891 for (int gt = 0; gt < ngtc; gt++)
893 real reft = std::max<real>(0, opts->ref_t[gt]);
894 if ((opts->tau_t[gt] > 0)
895 && (reft > 0)) /* tau_t or ref_t = 0 means that no randomization is done */
897 randomize_group[gt] = true;
898 boltzfac[gt] = gmx::c_boltz * opts->ref_t[gt];
902 randomize_group[gt] = false;
908 void Update::Impl::update_temperature_constants(const t_inputrec& inputRecord)
910 if (inputRecord.eI == IntegrationAlgorithm::BD)
912 if (inputRecord.bd_fric != 0)
914 for (int gt = 0; gt < inputRecord.opts.ngtc; gt++)
916 sd_.bd_rf[gt] = std::sqrt(2.0 * gmx::c_boltz * inputRecord.opts.ref_t[gt]
917 / (inputRecord.bd_fric * inputRecord.delta_t));
922 for (int gt = 0; gt < inputRecord.opts.ngtc; gt++)
924 sd_.bd_rf[gt] = std::sqrt(2.0 * gmx::c_boltz * inputRecord.opts.ref_t[gt]);
928 if (inputRecord.eI == IntegrationAlgorithm::SD1)
930 for (int gt = 0; gt < inputRecord.opts.ngtc; gt++)
932 real kT = gmx::c_boltz * inputRecord.opts.ref_t[gt];
933 /* The mass is accounted for later, since this differs per atom */
934 sd_.sdsig[gt].V = std::sqrt(kT * (1 - sd_.sdc[gt].em * sd_.sdc[gt].em));
939 Update::Impl::Impl(const t_inputrec& inputRecord, BoxDeformation* boxDeformation) :
941 deform_(boxDeformation)
943 update_temperature_constants(inputRecord);
944 xp_.resizeWithPadding(0);
947 void Update::setNumAtoms(int numAtoms)
950 impl_->xp()->resizeWithPadding(numAtoms);
953 /*! \brief Sets the SD update type */
954 enum class SDUpdate : int
957 FrictionAndNoiseOnly,
961 /*! \brief SD integrator update
963 * Two phases are required in the general case of a constrained
964 * update, the first phase from the contribution of forces, before
965 * applying constraints, and then a second phase applying the friction
966 * and noise, and then further constraining. For details, see
969 * Without constraints, the two phases can be combined, for
972 * Thus three instantiations of this templated function will be made,
973 * two with only one contribution, and one with both contributions. */
974 template<SDUpdate updateType>
975 static void doSDUpdateGeneral(const gmx_stochd_t& sd,
979 const ivec nFreeze[],
980 const real invmass[],
981 const ParticleType ptype[],
982 const unsigned short cFREEZE[],
983 const unsigned short cTC[],
992 // cTC and cFREEZE can be nullptr any time, but various
993 // instantiations do not make sense with particular pointer
995 if (updateType == SDUpdate::ForcesOnly)
997 GMX_ASSERT(f != nullptr, "SD update with only forces requires forces");
998 GMX_ASSERT(cTC == nullptr, "SD update with only forces cannot handle temperature groups");
1000 if (updateType == SDUpdate::FrictionAndNoiseOnly)
1002 GMX_ASSERT(f == nullptr, "SD update with only noise cannot handle forces");
1004 if (updateType == SDUpdate::Combined)
1006 GMX_ASSERT(f != nullptr, "SD update with forces and noise requires forces");
1009 // Even 0 bits internal counter gives 2x64 ints (more than enough for three table lookups)
1010 gmx::ThreeFry2x64<0> rng(seed, gmx::RandomDomain::UpdateCoordinates);
1011 gmx::TabulatedNormalDistribution<real, 14> dist;
1013 for (int n = start; n < nrend; n++)
1015 int globalAtomIndex = gatindex ? gatindex[n] : n;
1016 rng.restart(step, globalAtomIndex);
1019 real inverseMass = invmass[n];
1020 real invsqrtMass = std::sqrt(inverseMass);
1022 int freezeGroup = cFREEZE ? cFREEZE[n] : 0;
1023 int temperatureGroup = cTC ? cTC[n] : 0;
1025 for (int d = 0; d < DIM; d++)
1027 if ((ptype[n] != ParticleType::Shell) && !nFreeze[freezeGroup][d])
1029 if (updateType == SDUpdate::ForcesOnly)
1031 real vn = v[n][d] + inverseMass * f[n][d] * dt;
1033 // Simple position update.
1034 xprime[n][d] = x[n][d] + v[n][d] * dt;
1036 else if (updateType == SDUpdate::FrictionAndNoiseOnly)
1039 v[n][d] = (vn * sd.sdc[temperatureGroup].em
1040 + invsqrtMass * sd.sdsig[temperatureGroup].V * dist(rng));
1041 // The previous phase already updated the
1042 // positions with a full v*dt term that must
1043 // now be half removed.
1044 xprime[n][d] = xprime[n][d] + 0.5 * (v[n][d] - vn) * dt;
1048 real vn = v[n][d] + inverseMass * f[n][d] * dt;
1049 v[n][d] = (vn * sd.sdc[temperatureGroup].em
1050 + invsqrtMass * sd.sdsig[temperatureGroup].V * dist(rng));
1051 // Here we include half of the friction+noise
1052 // update of v into the position update.
1053 xprime[n][d] = x[n][d] + 0.5 * (vn + v[n][d]) * dt;
1058 // When using constraints, the update is split into
1059 // two phases, but we only need to zero the update of
1060 // virtual, shell or frozen particles in at most one
1062 if (updateType != SDUpdate::FrictionAndNoiseOnly)
1065 xprime[n][d] = x[n][d];
1072 static void do_update_sd(int start,
1076 const rvec* gmx_restrict x,
1077 rvec* gmx_restrict xprime,
1078 rvec* gmx_restrict v,
1079 const rvec* gmx_restrict f,
1080 const ivec nFreeze[],
1081 const real invmass[],
1082 const ParticleType ptype[],
1083 const unsigned short cFREEZE[],
1084 const unsigned short cTC[],
1086 const t_commrec* cr,
1087 const gmx_stochd_t& sd,
1088 bool haveConstraints)
1090 if (haveConstraints)
1092 // With constraints, the SD update is done in 2 parts
1093 doSDUpdateGeneral<SDUpdate::ForcesOnly>(
1094 sd, start, nrend, dt, nFreeze, invmass, ptype, cFREEZE, nullptr, x, xprime, v, f, step, seed, nullptr);
1098 doSDUpdateGeneral<SDUpdate::Combined>(sd,
1113 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1117 static void do_update_bd(int start,
1121 const rvec* gmx_restrict x,
1122 rvec* gmx_restrict xprime,
1123 rvec* gmx_restrict v,
1124 const rvec* gmx_restrict f,
1125 const ivec nFreeze[],
1126 const real invmass[],
1127 const ParticleType ptype[],
1128 const unsigned short cFREEZE[],
1129 const unsigned short cTC[],
1130 real friction_coefficient,
1133 const int* gatindex)
1135 /* note -- these appear to be full step velocities . . . */
1140 // Use 1 bit of internal counters to give us 2*2 64-bits values per stream
1141 // Each 64-bit value is enough for 4 normal distribution table numbers.
1142 gmx::ThreeFry2x64<0> rng(seed, gmx::RandomDomain::UpdateCoordinates);
1143 gmx::TabulatedNormalDistribution<real, 14> dist;
1145 if (friction_coefficient != 0)
1147 invfr = 1.0 / friction_coefficient;
1150 for (n = start; (n < nrend); n++)
1152 int ng = gatindex ? gatindex[n] : n;
1154 rng.restart(step, ng);
1165 for (d = 0; (d < DIM); d++)
1167 if ((ptype[n] != ParticleType::Shell) && !nFreeze[gf][d])
1169 if (friction_coefficient != 0)
1171 vn = invfr * f[n][d] + rf[gt] * dist(rng);
1175 /* NOTE: invmass = 2/(mass*friction_constant*dt) */
1176 vn = 0.5 * invmass[n] * f[n][d] * dt
1177 + std::sqrt(0.5 * invmass[n]) * rf[gt] * dist(rng);
1181 xprime[n][d] = x[n][d] + vn * dt;
1186 xprime[n][d] = x[n][d];
1192 extern void init_ekinstate(ekinstate_t* ekinstate, const t_inputrec* ir)
1194 ekinstate->ekin_n = ir->opts.ngtc;
1195 snew(ekinstate->ekinh, ekinstate->ekin_n);
1196 snew(ekinstate->ekinf, ekinstate->ekin_n);
1197 snew(ekinstate->ekinh_old, ekinstate->ekin_n);
1198 ekinstate->ekinscalef_nhc.resize(ekinstate->ekin_n);
1199 ekinstate->ekinscaleh_nhc.resize(ekinstate->ekin_n);
1200 ekinstate->vscale_nhc.resize(ekinstate->ekin_n);
1201 ekinstate->dekindl = 0;
1202 ekinstate->mvcos = 0;
1203 ekinstate->hasReadEkinState = false;
1206 void update_ekinstate(ekinstate_t* ekinstate, const gmx_ekindata_t* ekind)
1210 for (i = 0; i < ekinstate->ekin_n; i++)
1212 copy_mat(ekind->tcstat[i].ekinh, ekinstate->ekinh[i]);
1213 copy_mat(ekind->tcstat[i].ekinf, ekinstate->ekinf[i]);
1214 copy_mat(ekind->tcstat[i].ekinh_old, ekinstate->ekinh_old[i]);
1215 ekinstate->ekinscalef_nhc[i] = ekind->tcstat[i].ekinscalef_nhc;
1216 ekinstate->ekinscaleh_nhc[i] = ekind->tcstat[i].ekinscaleh_nhc;
1217 ekinstate->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
1220 copy_mat(ekind->ekin, ekinstate->ekin_total);
1221 ekinstate->dekindl = ekind->dekindl;
1222 ekinstate->mvcos = ekind->cosacc.mvcos;
1225 void restore_ekinstate_from_state(const t_commrec* cr, gmx_ekindata_t* ekind, const ekinstate_t* ekinstate)
1231 for (i = 0; i < ekinstate->ekin_n; i++)
1233 copy_mat(ekinstate->ekinh[i], ekind->tcstat[i].ekinh);
1234 copy_mat(ekinstate->ekinf[i], ekind->tcstat[i].ekinf);
1235 copy_mat(ekinstate->ekinh_old[i], ekind->tcstat[i].ekinh_old);
1236 ekind->tcstat[i].ekinscalef_nhc = ekinstate->ekinscalef_nhc[i];
1237 ekind->tcstat[i].ekinscaleh_nhc = ekinstate->ekinscaleh_nhc[i];
1238 ekind->tcstat[i].vscale_nhc = ekinstate->vscale_nhc[i];
1241 copy_mat(ekinstate->ekin_total, ekind->ekin);
1243 ekind->dekindl = ekinstate->dekindl;
1244 ekind->cosacc.mvcos = ekinstate->mvcos;
1245 n = ekinstate->ekin_n;
1250 gmx_bcast(sizeof(n), &n, cr->mpi_comm_mygroup);
1251 for (i = 0; i < n; i++)
1253 gmx_bcast(DIM * DIM * sizeof(ekind->tcstat[i].ekinh[0][0]),
1254 ekind->tcstat[i].ekinh[0],
1255 cr->mpi_comm_mygroup);
1256 gmx_bcast(DIM * DIM * sizeof(ekind->tcstat[i].ekinf[0][0]),
1257 ekind->tcstat[i].ekinf[0],
1258 cr->mpi_comm_mygroup);
1259 gmx_bcast(DIM * DIM * sizeof(ekind->tcstat[i].ekinh_old[0][0]),
1260 ekind->tcstat[i].ekinh_old[0],
1261 cr->mpi_comm_mygroup);
1263 gmx_bcast(sizeof(ekind->tcstat[i].ekinscalef_nhc),
1264 &(ekind->tcstat[i].ekinscalef_nhc),
1265 cr->mpi_comm_mygroup);
1266 gmx_bcast(sizeof(ekind->tcstat[i].ekinscaleh_nhc),
1267 &(ekind->tcstat[i].ekinscaleh_nhc),
1268 cr->mpi_comm_mygroup);
1269 gmx_bcast(sizeof(ekind->tcstat[i].vscale_nhc), &(ekind->tcstat[i].vscale_nhc), cr->mpi_comm_mygroup);
1271 gmx_bcast(DIM * DIM * sizeof(ekind->ekin[0][0]), ekind->ekin[0], cr->mpi_comm_mygroup);
1273 gmx_bcast(sizeof(ekind->dekindl), &ekind->dekindl, cr->mpi_comm_mygroup);
1274 gmx_bcast(sizeof(ekind->cosacc.mvcos), &ekind->cosacc.mvcos, cr->mpi_comm_mygroup);
1278 void getThreadAtomRange(int numThreads, int threadIndex, int numAtoms, int* startAtom, int* endAtom)
1280 #if GMX_HAVE_SIMD_UPDATE
1281 constexpr int blockSize = GMX_SIMD_REAL_WIDTH;
1283 constexpr int blockSize = 1;
1285 int numBlocks = (numAtoms + blockSize - 1) / blockSize;
1287 *startAtom = ((numBlocks * threadIndex) / numThreads) * blockSize;
1288 *endAtom = ((numBlocks * (threadIndex + 1)) / numThreads) * blockSize;
1289 if (threadIndex == numThreads - 1)
1291 *endAtom = numAtoms;
1295 void Update::Impl::update_sd_second_half(const t_inputrec& inputRecord,
1298 const t_mdatoms* md,
1300 const t_commrec* cr,
1302 gmx_wallcycle_t wcycle,
1303 gmx::Constraints* constr,
1311 if (inputRecord.eI == IntegrationAlgorithm::SD1)
1313 int homenr = md->homenr;
1315 /* Cast delta_t from double to real to make the integrators faster.
1316 * The only reason for having delta_t double is to get accurate values
1317 * for t=delta_t*step when step is larger than float precision.
1318 * For integration dt the accuracy of real suffices, since with
1319 * integral += dt*integrand the increment is nearly always (much) smaller
1320 * than the integral (and the integrand has real precision).
1322 real dt = inputRecord.delta_t;
1324 wallcycle_start(wcycle, ewcUPDATE);
1326 int nth = gmx_omp_nthreads_get(emntUpdate);
1328 #pragma omp parallel for num_threads(nth) schedule(static)
1329 for (int th = 0; th < nth; th++)
1333 int start_th, end_th;
1334 getThreadAtomRange(nth, th, homenr, &start_th, &end_th);
1336 doSDUpdateGeneral<SDUpdate::FrictionAndNoiseOnly>(
1341 inputRecord.opts.nFreeze,
1346 state->x.rvec_array(),
1348 state->v.rvec_array(),
1351 inputRecord.ld_seed,
1352 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1354 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1356 inc_nrnb(nrnb, eNR_UPDATE, homenr);
1357 wallcycle_stop(wcycle, ewcUPDATE);
1359 /* Constrain the coordinates upd->xp for half a time step */
1360 bool computeVirial = false;
1361 constr->apply(do_log,
1366 state->x.arrayRefWithPadding(),
1367 xp_.arrayRefWithPadding(),
1370 state->lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Bonded)],
1372 state->v.arrayRefWithPadding(),
1375 ConstraintVariable::Positions);
1379 void Update::Impl::finish_update(const t_inputrec& inputRecord,
1380 const t_mdatoms* md,
1382 gmx_wallcycle_t wcycle,
1383 const bool haveConstraints)
1385 /* NOTE: Currently we always integrate to a temporary buffer and
1386 * then copy the results back here.
1389 wallcycle_start_nocount(wcycle, ewcUPDATE);
1391 const int homenr = md->homenr;
1392 auto xp = makeConstArrayRef(xp_).subArray(0, homenr);
1393 auto x = makeArrayRef(state->x).subArray(0, homenr);
1395 if (md->havePartiallyFrozenAtoms && haveConstraints)
1397 /* We have atoms that are frozen along some, but not all dimensions,
1398 * then constraints will have moved them also along the frozen dimensions.
1399 * To freeze such degrees of freedom we do not copy them back here.
1401 const ivec* nFreeze = inputRecord.opts.nFreeze;
1403 for (int i = 0; i < homenr; i++)
1405 const int g = md->cFREEZE[i];
1407 for (int d = 0; d < DIM; d++)
1409 if (nFreeze[g][d] == 0)
1418 /* We have no frozen atoms or fully frozen atoms which have not
1419 * been moved by the update, so we can simply copy all coordinates.
1421 int gmx_unused nth = gmx_omp_nthreads_get(emntUpdate);
1422 #pragma omp parallel for num_threads(nth) schedule(static)
1423 for (int i = 0; i < homenr; i++)
1425 // Trivial statement, does not throw
1430 wallcycle_stop(wcycle, ewcUPDATE);
1433 void Update::Impl::update_coords(const t_inputrec& inputRecord,
1435 const t_mdatoms* md,
1437 const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
1438 const t_fcdata& fcdata,
1439 const gmx_ekindata_t* ekind,
1442 const t_commrec* cr,
1443 const bool haveConstraints)
1445 /* Running the velocity half does nothing except for velocity verlet */
1446 if ((updatePart == etrtVELOCITY1 || updatePart == etrtVELOCITY2) && !EI_VV(inputRecord.eI))
1448 gmx_incons("update_coords called for velocity without VV integrator");
1451 int homenr = md->homenr;
1453 /* Cast to real for faster code, no loss in precision (see comment above) */
1454 real dt = inputRecord.delta_t;
1456 /* We need to update the NMR restraint history when time averaging is used */
1457 if (state->flags & enumValueToBitMask(StateEntry::DisreRm3Tav))
1459 update_disres_history(*fcdata.disres, &state->hist);
1461 if (state->flags & enumValueToBitMask(StateEntry::OrireDtav))
1463 update_orires_history(*fcdata.orires, &state->hist);
1466 /* ############# START The update of velocities and positions ######### */
1467 int nth = gmx_omp_nthreads_get(emntUpdate);
1469 #pragma omp parallel for num_threads(nth) schedule(static)
1470 for (int th = 0; th < nth; th++)
1474 int start_th, end_th;
1475 getThreadAtomRange(nth, th, homenr, &start_th, &end_th);
1477 const rvec* x_rvec = state->x.rvec_array();
1478 rvec* xp_rvec = xp_.rvec_array();
1479 rvec* v_rvec = state->v.rvec_array();
1480 const rvec* f_rvec = as_rvec_array(f.unpaddedConstArrayRef().data());
1482 switch (inputRecord.eI)
1484 case (IntegrationAlgorithm::MD):
1485 do_update_md(start_th,
1495 inputRecord.nsttcouple,
1496 inputRecord.nstpcouple,
1500 state->nosehoover_vxi.data(),
1503 case (IntegrationAlgorithm::SD1):
1504 do_update_sd(start_th,
1512 inputRecord.opts.nFreeze,
1517 inputRecord.ld_seed,
1522 case (IntegrationAlgorithm::BD):
1523 do_update_bd(start_th,
1531 inputRecord.opts.nFreeze,
1536 inputRecord.bd_fric,
1538 inputRecord.ld_seed,
1539 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1541 case (IntegrationAlgorithm::VV):
1542 case (IntegrationAlgorithm::VVAK):
1544 gmx_bool bExtended = (inputRecord.etc == TemperatureCoupling::NoseHoover
1545 || inputRecord.epc == PressureCoupling::ParrinelloRahman
1546 || inputRecord.epc == PressureCoupling::Mttk);
1548 /* assuming barostat coupled to group 0 */
1549 real alpha = 1.0 + DIM / static_cast<real>(inputRecord.opts.nrdf[0]);
1554 do_update_vv_vel(start_th,
1557 inputRecord.opts.nFreeze,
1568 do_update_vv_pos(start_th,
1571 inputRecord.opts.nFreeze,
1583 default: gmx_fatal(FARGS, "Don't know how to update coordinates");
1586 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1590 void Update::Impl::update_for_constraint_virial(const t_inputrec& inputRecord,
1591 const t_mdatoms& md,
1592 const t_state& state,
1593 const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
1594 const gmx_ekindata_t& ekind)
1596 GMX_ASSERT(inputRecord.eI == IntegrationAlgorithm::MD || inputRecord.eI == IntegrationAlgorithm::SD1,
1597 "Only leap-frog is supported here");
1599 // Cast to real for faster code, no loss in precision
1600 const real dt = inputRecord.delta_t;
1602 const int nth = gmx_omp_nthreads_get(emntUpdate);
1604 #pragma omp parallel for num_threads(nth) schedule(static)
1605 for (int th = 0; th < nth; th++)
1609 int start_th, end_th;
1610 getThreadAtomRange(nth, th, md.homenr, &start_th, &end_th);
1612 const rvec* x_rvec = state.x.rvec_array();
1613 rvec* xp_rvec = xp_.rvec_array();
1614 const rvec* v_rvec = state.v.rvec_array();
1615 const rvec* f_rvec = as_rvec_array(f.unpaddedConstArrayRef().data());
1617 doUpdateMDDoNotUpdateVelocities(
1618 start_th, end_th, dt, x_rvec, xp_rvec, v_rvec, f_rvec, md, ekind);
1620 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR