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,
126 bool havePartiallyFrozenAtoms,
127 gmx::ArrayRef<const ParticleType> ptype,
128 gmx::ArrayRef<const real> invMass,
129 gmx::ArrayRef<const rvec> invMassPerDim,
131 const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
132 const t_fcdata& fcdata,
133 const gmx_ekindata_t* ekind,
137 bool haveConstraints);
139 void finish_update(const t_inputrec& inputRecord,
142 gmx_wallcycle* wcycle,
143 bool haveConstraints);
145 void update_sd_second_half(const t_inputrec& inputRecord,
149 gmx::ArrayRef<const ParticleType> ptype,
150 gmx::ArrayRef<const real> invMass,
154 gmx_wallcycle* wcycle,
155 gmx::Constraints* constr,
159 void update_for_constraint_virial(const t_inputrec& inputRecord,
161 bool havePartiallyFrozenAtoms,
162 gmx::ArrayRef<const real> invmass,
163 gmx::ArrayRef<const rvec> invMassPerDim,
164 const t_state& state,
165 const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
166 const gmx_ekindata_t& ekind);
168 void update_temperature_constants(const t_inputrec& inputRecord);
170 const std::vector<bool>& getAndersenRandomizeGroup() const { return sd_.randomize_group; }
172 const std::vector<real>& getBoltzmanFactor() const { return sd_.boltzfac; }
174 PaddedVector<RVec>* xp() { return &xp_; }
176 BoxDeformation* deform() const { return deform_; }
178 //! Group index for freezing
179 gmx::ArrayRef<const unsigned short> cFREEZE_;
180 //! Group index for temperature coupling
181 gmx::ArrayRef<const unsigned short> cTC_;
184 //! stochastic dynamics struct
186 //! xprime for constraint algorithms
187 PaddedVector<RVec> xp_;
188 //! Box deformation handler (or nullptr if inactive).
189 BoxDeformation* deform_ = nullptr;
192 Update::Update(const t_inputrec& inputRecord, BoxDeformation* boxDeformation) :
193 impl_(new Impl(inputRecord, boxDeformation)){};
197 const std::vector<bool>& Update::getAndersenRandomizeGroup() const
199 return impl_->getAndersenRandomizeGroup();
202 const std::vector<real>& Update::getBoltzmanFactor() const
204 return impl_->getBoltzmanFactor();
207 PaddedVector<RVec>* Update::xp()
212 BoxDeformation* Update::deform() const
214 return impl_->deform();
217 void Update::update_coords(const t_inputrec& inputRecord,
220 const bool havePartiallyFrozenAtoms,
221 gmx::ArrayRef<const ParticleType> ptype,
222 gmx::ArrayRef<const real> invMass,
223 gmx::ArrayRef<const rvec> invMassPerDim,
225 const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
226 const t_fcdata& fcdata,
227 const gmx_ekindata_t* ekind,
231 const bool haveConstraints)
233 return impl_->update_coords(inputRecord,
236 havePartiallyFrozenAtoms,
250 void Update::finish_update(const t_inputrec& inputRecord,
253 gmx_wallcycle* wcycle,
254 const bool haveConstraints)
256 return impl_->finish_update(inputRecord, md, state, wcycle, haveConstraints);
259 void Update::update_sd_second_half(const t_inputrec& inputRecord,
263 gmx::ArrayRef<const ParticleType> ptype,
264 gmx::ArrayRef<const real> invMass,
268 gmx_wallcycle* wcycle,
269 gmx::Constraints* constr,
273 return impl_->update_sd_second_half(
274 inputRecord, step, dvdlambda, homenr, ptype, invMass, state, cr, nrnb, wcycle, constr, do_log, do_ene);
277 void Update::update_for_constraint_virial(const t_inputrec& inputRecord,
279 const bool havePartiallyFrozenAtoms,
280 gmx::ArrayRef<const real> invmass,
281 gmx::ArrayRef<const rvec> invMassPerDim,
282 const t_state& state,
283 const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
284 const gmx_ekindata_t& ekind)
286 return impl_->update_for_constraint_virial(
287 inputRecord, homenr, havePartiallyFrozenAtoms, invmass, invMassPerDim, state, f, ekind);
290 void Update::update_temperature_constants(const t_inputrec& inputRecord)
292 return impl_->update_temperature_constants(inputRecord);
295 /*! \brief Sets whether we store the updated velocities */
296 enum class StoreUpdatedVelocities
298 yes, //!< Store the updated velocities
299 no //!< Do not store the updated velocities
302 /*! \brief Sets the number of different temperature coupling values */
303 enum class NumTempScaleValues
305 single, //!< Single T-scaling value (either one group or all values =1)
306 multiple //!< Multiple T-scaling values, need to use T-group indices
309 /*! \brief Sets if to apply no or diagonal Parrinello-Rahman pressure scaling
311 * Note that this enum is only used in updateMDLeapfrogSimple(), which does
312 * not handle fully anistropic Parrinello-Rahman scaling, so we only have
313 * options \p no and \p diagonal here and no anistropic option.
315 enum class ApplyParrinelloRahmanVScaling
317 no, //!< Do not apply velocity scaling (not a PR-coupling run or step)
318 diagonal //!< Apply velocity scaling using a diagonal matrix
321 /*! \brief Integrate using leap-frog with T-scaling and optionally diagonal Parrinello-Rahman p-coupling
323 * \tparam storeUpdatedVelocities Tells whether we should store the updated velocities
324 * \tparam numTempScaleValues The number of different T-couple values
325 * \tparam applyPRVScaling Apply Parrinello-Rahman velocity scaling
326 * \param[in] start Index of first atom to update
327 * \param[in] nrend Last atom to update: \p nrend - 1
328 * \param[in] dt The time step
329 * \param[in] dtPressureCouple Time step for pressure coupling
330 * \param[in] invMassPerDim 1/mass per atom and dimension
331 * \param[in] tcstat Temperature coupling information
332 * \param[in] cTC T-coupling group index per atom
333 * \param[in] pRVScaleMatrixDiagonal Parrinello-Rahman v-scale matrix diagonal
334 * \param[in] x Input coordinates
335 * \param[out] xprime Updated coordinates
336 * \param[inout] v Velocities, type either rvec* or const rvec*
337 * \param[in] f Forces
339 * We expect this template to get good SIMD acceleration by most compilers,
340 * unlike the more complex general template.
341 * Note that we might get even better SIMD acceleration when we introduce
342 * aligned (and padded) memory, possibly with some hints for the compilers.
344 template<StoreUpdatedVelocities storeUpdatedVelocities, NumTempScaleValues numTempScaleValues, ApplyParrinelloRahmanVScaling applyPRVScaling, typename VelocityType>
345 static std::enable_if_t<std::is_same<VelocityType, rvec>::value || std::is_same<VelocityType, const rvec>::value, void>
346 updateMDLeapfrogSimple(int start,
349 real dtPressureCouple,
350 gmx::ArrayRef<const rvec> invMassPerDim,
351 gmx::ArrayRef<const t_grp_tcstat> tcstat,
352 gmx::ArrayRef<const unsigned short> cTC,
353 const rvec pRVScaleMatrixDiagonal,
354 const rvec* gmx_restrict x,
355 rvec* gmx_restrict xprime,
356 VelocityType* gmx_restrict v,
357 const rvec* gmx_restrict f)
361 if (numTempScaleValues == NumTempScaleValues::single)
363 lambdaGroup = tcstat[0].lambda;
366 for (int a = start; a < nrend; a++)
368 if (numTempScaleValues == NumTempScaleValues::multiple)
370 lambdaGroup = tcstat[cTC[a]].lambda;
373 for (int d = 0; d < DIM; d++)
375 /* Note that using rvec invMassPerDim results in more efficient
376 * SIMD code, but this increases the cache pressure.
377 * For large systems with PME on the CPU this slows down the
378 * (then already slow) update by 20%. If all data remains in cache,
379 * using rvec is much faster.
381 real vNew = lambdaGroup * v[a][d] + f[a][d] * invMassPerDim[a][d] * dt;
383 if (applyPRVScaling == ApplyParrinelloRahmanVScaling::diagonal)
385 vNew -= dtPressureCouple * pRVScaleMatrixDiagonal[d] * v[a][d];
387 // TODO: Remove NOLINTs once clang-tidy is updated to v11, it should be able to handle constexpr.
388 if constexpr (storeUpdatedVelocities == StoreUpdatedVelocities::yes) // NOLINT // NOLINTNEXTLINE
392 xprime[a][d] = x[a][d] + vNew * dt; // NOLINT(readability-misleading-indentation)
397 #if GMX_SIMD && GMX_SIMD_HAVE_REAL
398 # define GMX_HAVE_SIMD_UPDATE 1
400 # define GMX_HAVE_SIMD_UPDATE 0
403 #if GMX_HAVE_SIMD_UPDATE
405 /*! \brief Load (aligned) the contents of GMX_SIMD_REAL_WIDTH rvec elements sequentially into 3 SIMD registers
407 * The loaded output is:
408 * \p r0: { r[index][XX], r[index][YY], ... }
410 * \p r2: { ..., r[index+GMX_SIMD_REAL_WIDTH-1][YY], r[index+GMX_SIMD_REAL_WIDTH-1][ZZ] }
412 * \param[in] r Real to an rvec array, has to be aligned to SIMD register width
413 * \param[in] index Index of the first rvec triplet of reals to load
414 * \param[out] r0 Pointer to first SIMD register
415 * \param[out] r1 Pointer to second SIMD register
416 * \param[out] r2 Pointer to third SIMD register
418 static inline void simdLoadRvecs(const rvec* r, int index, SimdReal* r0, SimdReal* r1, SimdReal* r2)
420 const real* realPtr = r[index];
422 GMX_ASSERT(isSimdAligned(realPtr), "Pointer should be SIMD aligned");
424 *r0 = simdLoad(realPtr + 0 * GMX_SIMD_REAL_WIDTH);
425 *r1 = simdLoad(realPtr + 1 * GMX_SIMD_REAL_WIDTH);
426 *r2 = simdLoad(realPtr + 2 * GMX_SIMD_REAL_WIDTH);
429 /*! \brief Store (aligned) 3 SIMD registers sequentially to GMX_SIMD_REAL_WIDTH rvec elements
431 * The stored output is:
432 * \p r[index] = { { r0[0], r0[1], ... }
434 * \p r[index+GMX_SIMD_REAL_WIDTH-1] = { ... , r2[GMX_SIMD_REAL_WIDTH-2], r2[GMX_SIMD_REAL_WIDTH-1] }
436 * \param[out] r Pointer to an rvec array, has to be aligned to SIMD register width
437 * \param[in] index Index of the first rvec triplet of reals to store to
438 * \param[in] r0 First SIMD register
439 * \param[in] r1 Second SIMD register
440 * \param[in] r2 Third SIMD register
442 static inline void simdStoreRvecs(rvec* r, int index, SimdReal r0, SimdReal r1, SimdReal r2)
444 real* realPtr = r[index];
446 GMX_ASSERT(isSimdAligned(realPtr), "Pointer should be SIMD aligned");
448 store(realPtr + 0 * GMX_SIMD_REAL_WIDTH, r0);
449 store(realPtr + 1 * GMX_SIMD_REAL_WIDTH, r1);
450 store(realPtr + 2 * GMX_SIMD_REAL_WIDTH, r2);
453 /*! \brief Integrate using leap-frog with single group T-scaling and SIMD
455 * \tparam storeUpdatedVelocities Tells whether we should store the updated velocities
456 * \param[in] start Index of first atom to update
457 * \param[in] nrend Last atom to update: \p nrend - 1
458 * \param[in] dt The time step
459 * \param[in] invMass 1/mass per atom
460 * \param[in] tcstat Temperature coupling information
461 * \param[in] x Input coordinates
462 * \param[out] xprime Updated coordinates
463 * \param[inout] v Velocities, type either rvec* or const rvec*
464 * \param[in] f Forces
466 template<StoreUpdatedVelocities storeUpdatedVelocities, typename VelocityType>
467 static std::enable_if_t<std::is_same<VelocityType, rvec>::value || std::is_same<VelocityType, const rvec>::value, void>
468 updateMDLeapfrogSimpleSimd(int start,
471 gmx::ArrayRef<const real> invMass,
472 gmx::ArrayRef<const t_grp_tcstat> tcstat,
473 const rvec* gmx_restrict x,
474 rvec* gmx_restrict xprime,
475 VelocityType* gmx_restrict v,
476 const rvec* gmx_restrict f)
478 SimdReal timestep(dt);
479 SimdReal lambdaSystem(tcstat[0].lambda);
481 /* We declare variables here, since code is often slower when declaring them inside the loop */
483 /* Note: We should implement a proper PaddedVector, so we don't need this check */
484 GMX_ASSERT(isSimdAligned(invMass.data()), "invMass should be aligned");
486 for (int a = start; a < nrend; a += GMX_SIMD_REAL_WIDTH)
488 SimdReal invMass0, invMass1, invMass2;
489 expandScalarsToTriplets(simdLoad(invMass.data() + a), &invMass0, &invMass1, &invMass2);
493 simdLoadRvecs(v, a, &v0, &v1, &v2);
494 simdLoadRvecs(f, a, &f0, &f1, &f2);
496 v0 = fma(f0 * invMass0, timestep, lambdaSystem * v0);
497 v1 = fma(f1 * invMass1, timestep, lambdaSystem * v1);
498 v2 = fma(f2 * invMass2, timestep, lambdaSystem * v2);
500 // TODO: Remove NOLINTs once clang-tidy is updated to v11, it should be able to handle constexpr.
501 if constexpr (storeUpdatedVelocities == StoreUpdatedVelocities::yes) // NOLINT // NOLINTNEXTLINE
503 simdStoreRvecs(v, a, v0, v1, v2);
506 SimdReal x0, x1, x2; // NOLINT(readability-misleading-indentation)
507 simdLoadRvecs(x, a, &x0, &x1, &x2);
509 SimdReal xprime0 = fma(v0, timestep, x0);
510 SimdReal xprime1 = fma(v1, timestep, x1);
511 SimdReal xprime2 = fma(v2, timestep, x2);
513 simdStoreRvecs(xprime, a, xprime0, xprime1, xprime2);
517 #endif // GMX_HAVE_SIMD_UPDATE
519 /*! \brief Sets the NEMD acceleration type */
520 enum class AccelerationType
526 /*! \brief Integrate using leap-frog with support for everything.
528 * \tparam accelerationType Type of NEMD acceleration.
529 * \param[in] start Index of first atom to update.
530 * \param[in] nrend Last atom to update: \p nrend - 1.
531 * \param[in] doNoseHoover If to apply Nose-Hoover T-coupling.
532 * \param[in] dt The time step.
533 * \param[in] dtPressureCouple Time step for pressure coupling, is 0 when no pressure
534 * coupling should be applied at this step.
535 * \param[in] cTC Temperature coupling group indices
536 * \param[in] invMassPerDim Inverse mass per dimension
537 * \param[in] ekind Kinetic energy data.
538 * \param[in] box The box dimensions.
539 * \param[in] x Input coordinates.
540 * \param[out] xprime Updated coordinates.
541 * \param[inout] v Velocities.
542 * \param[in] f Forces.
543 * \param[in] nh_vxi Nose-Hoover velocity scaling factors.
544 * \param[in] nsttcouple Frequency of the temperature coupling steps.
545 * \param[in] M Parrinello-Rahman scaling matrix.
547 template<AccelerationType accelerationType>
548 static void updateMDLeapfrogGeneral(int start,
552 real dtPressureCouple,
553 gmx::ArrayRef<const unsigned short> cTC,
554 gmx::ArrayRef<const rvec> invMassPerDim,
555 const gmx_ekindata_t* ekind,
557 const rvec* gmx_restrict x,
558 rvec* gmx_restrict xprime,
559 rvec* gmx_restrict v,
560 const rvec* gmx_restrict f,
561 const double* gmx_restrict nh_vxi,
562 const int nsttcouple,
565 /* This is a version of the leap-frog integrator that supports
566 * all combinations of T-coupling, P-coupling and NEMD.
567 * Nose-Hoover uses the reversible leap-frog integrator from
568 * Holian et al. Phys Rev E 52(3) : 2338, 1995
571 gmx::ArrayRef<const t_grp_tcstat> tcstat = ekind->tcstat;
573 /* Initialize group values, changed later when multiple groups are used */
576 real omega_Z = 2 * static_cast<real>(M_PI) / box[ZZ][ZZ];
578 for (int n = start; n < nrend; n++)
584 real lg = tcstat[gt].lambda;
587 real cosineZ, vCosine;
588 switch (accelerationType)
590 case AccelerationType::none: copy_rvec(v[n], vRel); break;
591 case AccelerationType::cosine:
592 cosineZ = std::cos(x[n][ZZ] * omega_Z);
593 vCosine = cosineZ * ekind->cosacc.vcos;
594 /* Avoid scaling the cosine profile velocity */
595 copy_rvec(v[n], vRel);
603 /* Here we account for multiple time stepping, by increasing
604 * the Nose-Hoover correction by nsttcouple
605 * TODO: This can be pre-computed.
607 factorNH = 0.5 * nsttcouple * dt * nh_vxi[gt];
610 for (int d = 0; d < DIM; d++)
612 real vNew = (lg * vRel[d]
613 + (f[n][d] * invMassPerDim[n][d] * dt - factorNH * vRel[d]
614 - dtPressureCouple * iprod(M[d], vRel)))
616 switch (accelerationType)
618 case AccelerationType::none: break;
619 case AccelerationType::cosine:
622 /* Add back the mean velocity and apply acceleration */
623 vNew += vCosine + cosineZ * ekind->cosacc.cos_accel * dt;
628 xprime[n][d] = x[n][d] + vNew * dt;
633 /*! \brief Handles the Leap-frog MD x and v integration */
634 static void do_update_md(int start,
638 const rvec* gmx_restrict x,
639 rvec* gmx_restrict xprime,
640 rvec* gmx_restrict v,
641 const rvec* gmx_restrict f,
642 const TemperatureCoupling etc,
643 const PressureCoupling epc,
644 const int nsttcouple,
645 const int nstpcouple,
646 gmx::ArrayRef<const unsigned short> cTC,
647 gmx::ArrayRef<const real> gmx_unused invmass,
648 gmx::ArrayRef<const rvec> invMassPerDim,
649 const gmx_ekindata_t* ekind,
651 const double* gmx_restrict nh_vxi,
653 bool gmx_unused havePartiallyFrozenAtoms)
655 GMX_ASSERT(nrend == start || xprime != x,
656 "For SIMD optimization certain compilers need to have xprime != x");
658 /* Note: Berendsen pressure scaling is handled after do_update_md() */
660 (etc != TemperatureCoupling::No && do_per_step(step + nsttcouple - 1, nsttcouple));
661 bool doNoseHoover = (etc == TemperatureCoupling::NoseHoover && doTempCouple);
662 bool doParrinelloRahman = (epc == PressureCoupling::ParrinelloRahman
663 && do_per_step(step + nstpcouple - 1, nstpcouple));
664 bool doPROffDiagonal = (doParrinelloRahman && (M[YY][XX] != 0 || M[ZZ][XX] != 0 || M[ZZ][YY] != 0));
666 real dtPressureCouple = (doParrinelloRahman ? nstpcouple * dt : 0);
668 /* NEMD cosine acceleration is applied in updateMDLeapFrogGeneral */
669 bool doAcceleration = (ekind->cosacc.cos_accel != 0);
671 if (doNoseHoover || doPROffDiagonal || doAcceleration)
674 if (!doParrinelloRahman)
676 /* We should not apply PR scaling at this step */
686 updateMDLeapfrogGeneral<AccelerationType::none>(start,
705 updateMDLeapfrogGeneral<AccelerationType::cosine>(start,
725 /* We determine if we have a single T-coupling lambda value for all
726 * atoms. That allows for better SIMD acceleration in the template.
727 * If we do not do temperature coupling (in the run or this step),
728 * all scaling values are 1, so we effectively have a single value.
730 bool haveSingleTempScaleValue = (!doTempCouple || ekind->ngtc == 1);
732 /* Extract some pointers needed by all cases */
733 gmx::ArrayRef<const t_grp_tcstat> tcstat = ekind->tcstat;
735 if (doParrinelloRahman)
737 GMX_ASSERT(!doPROffDiagonal,
738 "updateMDLeapfrogSimple only support diagonal Parrinello-Rahman scaling "
742 for (int d = 0; d < DIM; d++)
747 if (haveSingleTempScaleValue)
749 updateMDLeapfrogSimple<StoreUpdatedVelocities::yes, NumTempScaleValues::single, ApplyParrinelloRahmanVScaling::diagonal>(
750 start, nrend, dt, dtPressureCouple, invMassPerDim, tcstat, cTC, diagM, x, xprime, v, f);
754 updateMDLeapfrogSimple<StoreUpdatedVelocities::yes, NumTempScaleValues::multiple, ApplyParrinelloRahmanVScaling::diagonal>(
755 start, nrend, dt, dtPressureCouple, invMassPerDim, tcstat, cTC, diagM, x, xprime, v, f);
760 if (haveSingleTempScaleValue)
762 /* Note that modern compilers are pretty good at vectorizing
763 * updateMDLeapfrogSimple(). But the SIMD version will still
764 * be faster because invMass lowers the cache pressure
765 * compared to invMassPerDim.
767 #if GMX_HAVE_SIMD_UPDATE
768 /* Check if we can use invmass instead of invMassPerDim */
769 if (!havePartiallyFrozenAtoms)
771 updateMDLeapfrogSimpleSimd<StoreUpdatedVelocities::yes>(
772 start, nrend, dt, invmass, tcstat, x, xprime, v, f);
777 updateMDLeapfrogSimple<StoreUpdatedVelocities::yes, NumTempScaleValues::single, ApplyParrinelloRahmanVScaling::no>(
778 start, nrend, dt, dtPressureCouple, invMassPerDim, tcstat, cTC, nullptr, x, xprime, v, f);
783 updateMDLeapfrogSimple<StoreUpdatedVelocities::yes, NumTempScaleValues::multiple, ApplyParrinelloRahmanVScaling::no>(
784 start, nrend, dt, dtPressureCouple, invMassPerDim, tcstat, cTC, nullptr, x, xprime, v, f);
789 /*! \brief Handles the Leap-frog MD x and v integration */
790 static void doUpdateMDDoNotUpdateVelocities(int start,
793 const rvec* gmx_restrict x,
794 rvec* gmx_restrict xprime,
795 const rvec* gmx_restrict v,
796 const rvec* gmx_restrict f,
797 bool gmx_unused havePartiallyFrozenAtoms,
798 gmx::ArrayRef<const real> gmx_unused invmass,
799 gmx::ArrayRef<const rvec> invMassPerDim,
800 const gmx_ekindata_t& ekind)
802 GMX_ASSERT(nrend == start || xprime != x,
803 "For SIMD optimization certain compilers need to have xprime != x");
805 gmx::ArrayRef<const t_grp_tcstat> tcstat = ekind.tcstat;
807 /* Check if we can use invmass instead of invMassPerDim */
808 #if GMX_HAVE_SIMD_UPDATE
809 if (!havePartiallyFrozenAtoms)
811 updateMDLeapfrogSimpleSimd<StoreUpdatedVelocities::no>(
812 start, nrend, dt, invmass, tcstat, x, xprime, v, f);
817 updateMDLeapfrogSimple<StoreUpdatedVelocities::no, NumTempScaleValues::single, ApplyParrinelloRahmanVScaling::no>(
818 start, nrend, dt, dt, invMassPerDim, tcstat, gmx::ArrayRef<const unsigned short>(), nullptr, x, xprime, v, f);
822 static void do_update_vv_vel(int start,
825 gmx::ArrayRef<const ivec> nFreeze,
826 gmx::ArrayRef<const real> invmass,
827 gmx::ArrayRef<const ParticleType> ptype,
828 gmx::ArrayRef<const unsigned short> cFREEZE,
841 g = 0.25 * dt * veta * alpha;
843 mv2 = gmx::series_sinhx(g);
850 for (n = start; n < nrend; n++)
852 real w_dt = invmass[n] * dt;
853 if (!cFREEZE.empty())
858 for (d = 0; d < DIM; d++)
860 if ((ptype[n] != ParticleType::Shell) && !nFreeze[gf][d])
862 v[n][d] = mv1 * (mv1 * v[n][d] + 0.5 * (w_dt * mv2 * f[n][d]));
870 } /* do_update_vv_vel */
872 static void do_update_vv_pos(int start,
875 gmx::ArrayRef<const ivec> nFreeze,
876 gmx::ArrayRef<const ParticleType> ptype,
877 gmx::ArrayRef<const unsigned short> cFREEZE,
888 /* Would it make more sense if Parrinello-Rahman was put here? */
893 mr2 = gmx::series_sinhx(g);
901 for (n = start; n < nrend; n++)
904 if (!cFREEZE.empty())
909 for (d = 0; d < DIM; d++)
911 if ((ptype[n] != ParticleType::Shell) && !nFreeze[gf][d])
913 xprime[n][d] = mr1 * (mr1 * x[n][d] + mr2 * dt * v[n][d]);
917 xprime[n][d] = x[n][d];
921 } /* do_update_vv_pos */
923 gmx_stochd_t::gmx_stochd_t(const t_inputrec& inputRecord)
925 const t_grpopts* opts = &inputRecord.opts;
926 const int ngtc = opts->ngtc;
928 if (inputRecord.eI == IntegrationAlgorithm::BD)
932 else if (EI_SD(inputRecord.eI))
937 for (int gt = 0; gt < ngtc; gt++)
939 if (opts->tau_t[gt] > 0)
941 sdc[gt].em = std::exp(-inputRecord.delta_t / opts->tau_t[gt]);
945 /* No friction and noise on this group */
950 else if (ETC_ANDERSEN(inputRecord.etc))
952 randomize_group.resize(ngtc);
953 boltzfac.resize(ngtc);
955 /* for now, assume that all groups, if randomized, are randomized at the same rate, i.e. tau_t is the same. */
956 /* since constraint groups don't necessarily match up with temperature groups! This is checked in readir.c */
958 for (int gt = 0; gt < ngtc; gt++)
960 real reft = std::max<real>(0, opts->ref_t[gt]);
961 if ((opts->tau_t[gt] > 0)
962 && (reft > 0)) /* tau_t or ref_t = 0 means that no randomization is done */
964 randomize_group[gt] = true;
965 boltzfac[gt] = gmx::c_boltz * opts->ref_t[gt];
969 randomize_group[gt] = false;
975 void Update::Impl::update_temperature_constants(const t_inputrec& inputRecord)
977 if (inputRecord.eI == IntegrationAlgorithm::BD)
979 if (inputRecord.bd_fric != 0)
981 for (int gt = 0; gt < inputRecord.opts.ngtc; gt++)
983 sd_.bd_rf[gt] = std::sqrt(2.0 * gmx::c_boltz * inputRecord.opts.ref_t[gt]
984 / (inputRecord.bd_fric * inputRecord.delta_t));
989 for (int gt = 0; gt < inputRecord.opts.ngtc; gt++)
991 sd_.bd_rf[gt] = std::sqrt(2.0 * gmx::c_boltz * inputRecord.opts.ref_t[gt]);
995 if (inputRecord.eI == IntegrationAlgorithm::SD1)
997 for (int gt = 0; gt < inputRecord.opts.ngtc; gt++)
999 real kT = gmx::c_boltz * inputRecord.opts.ref_t[gt];
1000 /* The mass is accounted for later, since this differs per atom */
1001 sd_.sdsig[gt].V = std::sqrt(kT * (1 - sd_.sdc[gt].em * sd_.sdc[gt].em));
1006 Update::Impl::Impl(const t_inputrec& inputRecord, BoxDeformation* boxDeformation) :
1008 deform_(boxDeformation)
1010 update_temperature_constants(inputRecord);
1011 xp_.resizeWithPadding(0);
1014 void Update::updateAfterPartition(int numAtoms,
1015 gmx::ArrayRef<const unsigned short> cFREEZE,
1016 gmx::ArrayRef<const unsigned short> cTC)
1019 impl_->xp()->resizeWithPadding(numAtoms);
1020 impl_->cFREEZE_ = cFREEZE;
1024 /*! \brief Sets the SD update type */
1025 enum class SDUpdate : int
1028 FrictionAndNoiseOnly,
1032 /*! \brief SD integrator update
1034 * Two phases are required in the general case of a constrained
1035 * update, the first phase from the contribution of forces, before
1036 * applying constraints, and then a second phase applying the friction
1037 * and noise, and then further constraining. For details, see
1040 * Without constraints, the two phases can be combined, for
1043 * Thus three instantiations of this templated function will be made,
1044 * two with only one contribution, and one with both contributions. */
1045 template<SDUpdate updateType>
1046 static void doSDUpdateGeneral(const gmx_stochd_t& sd,
1050 gmx::ArrayRef<const ivec> nFreeze,
1051 gmx::ArrayRef<const real> invmass,
1052 gmx::ArrayRef<const ParticleType> ptype,
1053 gmx::ArrayRef<const unsigned short> cFREEZE,
1054 gmx::ArrayRef<const unsigned short> cTC,
1061 const int* gatindex)
1063 // cTC and cFREEZE can be nullptr any time, but various
1064 // instantiations do not make sense with particular pointer
1066 if (updateType == SDUpdate::ForcesOnly)
1068 GMX_ASSERT(f != nullptr, "SD update with only forces requires forces");
1069 GMX_ASSERT(cTC.empty(), "SD update with only forces cannot handle temperature groups");
1071 if (updateType == SDUpdate::FrictionAndNoiseOnly)
1073 GMX_ASSERT(f == nullptr, "SD update with only noise cannot handle forces");
1075 if (updateType == SDUpdate::Combined)
1077 GMX_ASSERT(f != nullptr, "SD update with forces and noise requires forces");
1080 // Even 0 bits internal counter gives 2x64 ints (more than enough for three table lookups)
1081 gmx::ThreeFry2x64<0> rng(seed, gmx::RandomDomain::UpdateCoordinates);
1082 gmx::TabulatedNormalDistribution<real, 14> dist;
1084 for (int n = start; n < nrend; n++)
1086 int globalAtomIndex = gatindex ? gatindex[n] : n;
1087 rng.restart(step, globalAtomIndex);
1090 real inverseMass = invmass[n];
1091 real invsqrtMass = std::sqrt(inverseMass);
1093 int freezeGroup = !cFREEZE.empty() ? cFREEZE[n] : 0;
1094 int temperatureGroup = !cTC.empty() ? cTC[n] : 0;
1096 for (int d = 0; d < DIM; d++)
1098 if ((ptype[n] != ParticleType::Shell) && !nFreeze[freezeGroup][d])
1100 if (updateType == SDUpdate::ForcesOnly)
1102 real vn = v[n][d] + inverseMass * f[n][d] * dt;
1104 // Simple position update.
1105 xprime[n][d] = x[n][d] + v[n][d] * dt;
1107 else if (updateType == SDUpdate::FrictionAndNoiseOnly)
1110 v[n][d] = (vn * sd.sdc[temperatureGroup].em
1111 + invsqrtMass * sd.sdsig[temperatureGroup].V * dist(rng));
1112 // The previous phase already updated the
1113 // positions with a full v*dt term that must
1114 // now be half removed.
1115 xprime[n][d] = xprime[n][d] + 0.5 * (v[n][d] - vn) * dt;
1119 real vn = v[n][d] + inverseMass * f[n][d] * dt;
1120 v[n][d] = (vn * sd.sdc[temperatureGroup].em
1121 + invsqrtMass * sd.sdsig[temperatureGroup].V * dist(rng));
1122 // Here we include half of the friction+noise
1123 // update of v into the position update.
1124 xprime[n][d] = x[n][d] + 0.5 * (vn + v[n][d]) * dt;
1129 // When using constraints, the update is split into
1130 // two phases, but we only need to zero the update of
1131 // virtual, shell or frozen particles in at most one
1133 if (updateType != SDUpdate::FrictionAndNoiseOnly)
1136 xprime[n][d] = x[n][d];
1143 static void do_update_sd(int start,
1147 const rvec* gmx_restrict x,
1148 rvec* gmx_restrict xprime,
1149 rvec* gmx_restrict v,
1150 const rvec* gmx_restrict f,
1151 gmx::ArrayRef<const ivec> nFreeze,
1152 gmx::ArrayRef<const real> invmass,
1153 gmx::ArrayRef<const ParticleType> ptype,
1154 gmx::ArrayRef<const unsigned short> cFREEZE,
1155 gmx::ArrayRef<const unsigned short> cTC,
1157 const t_commrec* cr,
1158 const gmx_stochd_t& sd,
1159 bool haveConstraints)
1161 if (haveConstraints)
1163 // With constraints, the SD update is done in 2 parts
1164 doSDUpdateGeneral<SDUpdate::ForcesOnly>(sd,
1172 gmx::ArrayRef<const unsigned short>(),
1183 doSDUpdateGeneral<SDUpdate::Combined>(sd,
1198 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1202 static void do_update_bd(int start,
1206 const rvec* gmx_restrict x,
1207 rvec* gmx_restrict xprime,
1208 rvec* gmx_restrict v,
1209 const rvec* gmx_restrict f,
1210 gmx::ArrayRef<const ivec> nFreeze,
1211 gmx::ArrayRef<const real> invmass,
1212 gmx::ArrayRef<const ParticleType> ptype,
1213 gmx::ArrayRef<const unsigned short> cFREEZE,
1214 gmx::ArrayRef<const unsigned short> cTC,
1215 real friction_coefficient,
1218 const int* gatindex)
1220 /* note -- these appear to be full step velocities . . . */
1225 // Use 1 bit of internal counters to give us 2*2 64-bits values per stream
1226 // Each 64-bit value is enough for 4 normal distribution table numbers.
1227 gmx::ThreeFry2x64<0> rng(seed, gmx::RandomDomain::UpdateCoordinates);
1228 gmx::TabulatedNormalDistribution<real, 14> dist;
1230 if (friction_coefficient != 0)
1232 invfr = 1.0 / friction_coefficient;
1235 for (n = start; (n < nrend); n++)
1237 int ng = gatindex ? gatindex[n] : n;
1239 rng.restart(step, ng);
1242 if (!cFREEZE.empty())
1250 for (d = 0; (d < DIM); d++)
1252 if ((ptype[n] != ParticleType::Shell) && !nFreeze[gf][d])
1254 if (friction_coefficient != 0)
1256 vn = invfr * f[n][d] + rf[gt] * dist(rng);
1260 /* NOTE: invmass = 2/(mass*friction_constant*dt) */
1261 vn = 0.5 * invmass[n] * f[n][d] * dt
1262 + std::sqrt(0.5 * invmass[n]) * rf[gt] * dist(rng);
1266 xprime[n][d] = x[n][d] + vn * dt;
1271 xprime[n][d] = x[n][d];
1277 extern void init_ekinstate(ekinstate_t* ekinstate, const t_inputrec* ir)
1279 ekinstate->ekin_n = ir->opts.ngtc;
1280 snew(ekinstate->ekinh, ekinstate->ekin_n);
1281 snew(ekinstate->ekinf, ekinstate->ekin_n);
1282 snew(ekinstate->ekinh_old, ekinstate->ekin_n);
1283 ekinstate->ekinscalef_nhc.resize(ekinstate->ekin_n);
1284 ekinstate->ekinscaleh_nhc.resize(ekinstate->ekin_n);
1285 ekinstate->vscale_nhc.resize(ekinstate->ekin_n);
1286 ekinstate->dekindl = 0;
1287 ekinstate->mvcos = 0;
1288 ekinstate->hasReadEkinState = false;
1291 void update_ekinstate(ekinstate_t* ekinstate, const gmx_ekindata_t* ekind)
1295 for (i = 0; i < ekinstate->ekin_n; i++)
1297 copy_mat(ekind->tcstat[i].ekinh, ekinstate->ekinh[i]);
1298 copy_mat(ekind->tcstat[i].ekinf, ekinstate->ekinf[i]);
1299 copy_mat(ekind->tcstat[i].ekinh_old, ekinstate->ekinh_old[i]);
1300 ekinstate->ekinscalef_nhc[i] = ekind->tcstat[i].ekinscalef_nhc;
1301 ekinstate->ekinscaleh_nhc[i] = ekind->tcstat[i].ekinscaleh_nhc;
1302 ekinstate->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
1305 copy_mat(ekind->ekin, ekinstate->ekin_total);
1306 ekinstate->dekindl = ekind->dekindl;
1307 ekinstate->mvcos = ekind->cosacc.mvcos;
1310 void restore_ekinstate_from_state(const t_commrec* cr, gmx_ekindata_t* ekind, const ekinstate_t* ekinstate)
1316 for (i = 0; i < ekinstate->ekin_n; i++)
1318 copy_mat(ekinstate->ekinh[i], ekind->tcstat[i].ekinh);
1319 copy_mat(ekinstate->ekinf[i], ekind->tcstat[i].ekinf);
1320 copy_mat(ekinstate->ekinh_old[i], ekind->tcstat[i].ekinh_old);
1321 ekind->tcstat[i].ekinscalef_nhc = ekinstate->ekinscalef_nhc[i];
1322 ekind->tcstat[i].ekinscaleh_nhc = ekinstate->ekinscaleh_nhc[i];
1323 ekind->tcstat[i].vscale_nhc = ekinstate->vscale_nhc[i];
1326 copy_mat(ekinstate->ekin_total, ekind->ekin);
1328 ekind->dekindl = ekinstate->dekindl;
1329 ekind->cosacc.mvcos = ekinstate->mvcos;
1330 n = ekinstate->ekin_n;
1335 gmx_bcast(sizeof(n), &n, cr->mpi_comm_mygroup);
1336 for (i = 0; i < n; i++)
1338 gmx_bcast(DIM * DIM * sizeof(ekind->tcstat[i].ekinh[0][0]),
1339 ekind->tcstat[i].ekinh[0],
1340 cr->mpi_comm_mygroup);
1341 gmx_bcast(DIM * DIM * sizeof(ekind->tcstat[i].ekinf[0][0]),
1342 ekind->tcstat[i].ekinf[0],
1343 cr->mpi_comm_mygroup);
1344 gmx_bcast(DIM * DIM * sizeof(ekind->tcstat[i].ekinh_old[0][0]),
1345 ekind->tcstat[i].ekinh_old[0],
1346 cr->mpi_comm_mygroup);
1348 gmx_bcast(sizeof(ekind->tcstat[i].ekinscalef_nhc),
1349 &(ekind->tcstat[i].ekinscalef_nhc),
1350 cr->mpi_comm_mygroup);
1351 gmx_bcast(sizeof(ekind->tcstat[i].ekinscaleh_nhc),
1352 &(ekind->tcstat[i].ekinscaleh_nhc),
1353 cr->mpi_comm_mygroup);
1354 gmx_bcast(sizeof(ekind->tcstat[i].vscale_nhc), &(ekind->tcstat[i].vscale_nhc), cr->mpi_comm_mygroup);
1356 gmx_bcast(DIM * DIM * sizeof(ekind->ekin[0][0]), ekind->ekin[0], cr->mpi_comm_mygroup);
1358 gmx_bcast(sizeof(ekind->dekindl), &ekind->dekindl, cr->mpi_comm_mygroup);
1359 gmx_bcast(sizeof(ekind->cosacc.mvcos), &ekind->cosacc.mvcos, cr->mpi_comm_mygroup);
1363 void getThreadAtomRange(int numThreads, int threadIndex, int numAtoms, int* startAtom, int* endAtom)
1365 #if GMX_HAVE_SIMD_UPDATE
1366 constexpr int blockSize = GMX_SIMD_REAL_WIDTH;
1368 constexpr int blockSize = 1;
1370 int numBlocks = (numAtoms + blockSize - 1) / blockSize;
1372 *startAtom = ((numBlocks * threadIndex) / numThreads) * blockSize;
1373 *endAtom = ((numBlocks * (threadIndex + 1)) / numThreads) * blockSize;
1374 if (threadIndex == numThreads - 1)
1376 *endAtom = numAtoms;
1380 void Update::Impl::update_sd_second_half(const t_inputrec& inputRecord,
1384 gmx::ArrayRef<const ParticleType> ptype,
1385 gmx::ArrayRef<const real> invMass,
1387 const t_commrec* cr,
1389 gmx_wallcycle* wcycle,
1390 gmx::Constraints* constr,
1398 if (inputRecord.eI == IntegrationAlgorithm::SD1)
1401 /* Cast delta_t from double to real to make the integrators faster.
1402 * The only reason for having delta_t double is to get accurate values
1403 * for t=delta_t*step when step is larger than float precision.
1404 * For integration dt the accuracy of real suffices, since with
1405 * integral += dt*integrand the increment is nearly always (much) smaller
1406 * than the integral (and the integrand has real precision).
1408 real dt = inputRecord.delta_t;
1410 wallcycle_start(wcycle, WallCycleCounter::Update);
1412 int nth = gmx_omp_nthreads_get(ModuleMultiThread::Update);
1414 #pragma omp parallel for num_threads(nth) schedule(static)
1415 for (int th = 0; th < nth; th++)
1419 int start_th, end_th;
1420 getThreadAtomRange(nth, th, homenr, &start_th, &end_th);
1422 doSDUpdateGeneral<SDUpdate::FrictionAndNoiseOnly>(
1427 gmx::arrayRefFromArray(inputRecord.opts.nFreeze, inputRecord.opts.ngfrz),
1432 state->x.rvec_array(),
1434 state->v.rvec_array(),
1437 inputRecord.ld_seed,
1438 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1440 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1442 inc_nrnb(nrnb, eNR_UPDATE, homenr);
1443 wallcycle_stop(wcycle, WallCycleCounter::Update);
1445 /* Constrain the coordinates upd->xp for half a time step */
1446 bool computeVirial = false;
1447 constr->apply(do_log,
1452 state->x.arrayRefWithPadding(),
1453 xp_.arrayRefWithPadding(),
1456 state->lambda[static_cast<int>(FreeEnergyPerturbationCouplingType::Bonded)],
1458 state->v.arrayRefWithPadding(),
1461 ConstraintVariable::Positions);
1465 void Update::Impl::finish_update(const t_inputrec& inputRecord,
1466 const t_mdatoms* md,
1468 gmx_wallcycle* wcycle,
1469 const bool haveConstraints)
1471 /* NOTE: Currently we always integrate to a temporary buffer and
1472 * then copy the results back here.
1475 wallcycle_start_nocount(wcycle, WallCycleCounter::Update);
1477 const int homenr = md->homenr;
1478 auto xp = makeConstArrayRef(xp_).subArray(0, homenr);
1479 auto x = makeArrayRef(state->x).subArray(0, homenr);
1481 if (md->havePartiallyFrozenAtoms && haveConstraints)
1483 /* We have atoms that are frozen along some, but not all dimensions,
1484 * then constraints will have moved them also along the frozen dimensions.
1485 * To freeze such degrees of freedom we do not copy them back here.
1487 const ivec* nFreeze = inputRecord.opts.nFreeze;
1489 for (int i = 0; i < homenr; i++)
1491 const int g = md->cFREEZE[i];
1493 for (int d = 0; d < DIM; d++)
1495 if (nFreeze[g][d] == 0)
1504 /* We have no frozen atoms or fully frozen atoms which have not
1505 * been moved by the update, so we can simply copy all coordinates.
1507 int gmx_unused nth = gmx_omp_nthreads_get(ModuleMultiThread::Update);
1508 #pragma omp parallel for num_threads(nth) schedule(static)
1509 for (int i = 0; i < homenr; i++)
1511 // Trivial statement, does not throw
1516 wallcycle_stop(wcycle, WallCycleCounter::Update);
1519 void Update::Impl::update_coords(const t_inputrec& inputRecord,
1522 bool havePartiallyFrozenAtoms,
1523 gmx::ArrayRef<const ParticleType> ptype,
1524 gmx::ArrayRef<const real> invMass,
1525 gmx::ArrayRef<const rvec> invMassPerDim,
1527 const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
1528 const t_fcdata& fcdata,
1529 const gmx_ekindata_t* ekind,
1532 const t_commrec* cr,
1533 const bool haveConstraints)
1535 /* Running the velocity half does nothing except for velocity verlet */
1536 if ((updatePart == etrtVELOCITY1 || updatePart == etrtVELOCITY2) && !EI_VV(inputRecord.eI))
1538 gmx_incons("update_coords called for velocity without VV integrator");
1541 /* Cast to real for faster code, no loss in precision (see comment above) */
1542 real dt = inputRecord.delta_t;
1544 /* We need to update the NMR restraint history when time averaging is used */
1545 if (state->flags & enumValueToBitMask(StateEntry::DisreRm3Tav))
1547 update_disres_history(*fcdata.disres, &state->hist);
1549 if (state->flags & enumValueToBitMask(StateEntry::OrireDtav))
1551 update_orires_history(*fcdata.orires, &state->hist);
1554 /* ############# START The update of velocities and positions ######### */
1555 int nth = gmx_omp_nthreads_get(ModuleMultiThread::Update);
1557 #pragma omp parallel for num_threads(nth) schedule(static)
1558 for (int th = 0; th < nth; th++)
1562 int start_th, end_th;
1563 getThreadAtomRange(nth, th, homenr, &start_th, &end_th);
1565 const rvec* x_rvec = state->x.rvec_array();
1566 rvec* xp_rvec = xp_.rvec_array();
1567 rvec* v_rvec = state->v.rvec_array();
1568 const rvec* f_rvec = as_rvec_array(f.unpaddedConstArrayRef().data());
1570 switch (inputRecord.eI)
1572 case (IntegrationAlgorithm::MD):
1573 do_update_md(start_th,
1583 inputRecord.nsttcouple,
1584 inputRecord.nstpcouple,
1590 state->nosehoover_vxi.data(),
1592 havePartiallyFrozenAtoms);
1594 case (IntegrationAlgorithm::SD1):
1595 do_update_sd(start_th,
1603 gmx::arrayRefFromArray(inputRecord.opts.nFreeze, inputRecord.opts.ngfrz),
1608 inputRecord.ld_seed,
1613 case (IntegrationAlgorithm::BD):
1614 do_update_bd(start_th,
1622 gmx::arrayRefFromArray(inputRecord.opts.nFreeze, inputRecord.opts.ngfrz),
1627 inputRecord.bd_fric,
1629 inputRecord.ld_seed,
1630 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1632 case (IntegrationAlgorithm::VV):
1633 case (IntegrationAlgorithm::VVAK):
1635 gmx_bool bExtended = (inputRecord.etc == TemperatureCoupling::NoseHoover
1636 || inputRecord.epc == PressureCoupling::ParrinelloRahman
1637 || inputRecord.epc == PressureCoupling::Mttk);
1639 /* assuming barostat coupled to group 0 */
1640 real alpha = 1.0 + DIM / static_cast<real>(inputRecord.opts.nrdf[0]);
1645 do_update_vv_vel(start_th,
1648 gmx::arrayRefFromArray(inputRecord.opts.nFreeze,
1649 inputRecord.opts.ngfrz),
1660 do_update_vv_pos(start_th,
1663 gmx::arrayRefFromArray(inputRecord.opts.nFreeze,
1664 inputRecord.opts.ngfrz),
1676 default: gmx_fatal(FARGS, "Don't know how to update coordinates");
1679 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1683 void Update::Impl::update_for_constraint_virial(const t_inputrec& inputRecord,
1685 bool havePartiallyFrozenAtoms,
1686 gmx::ArrayRef<const real> invmass,
1687 gmx::ArrayRef<const rvec> invMassPerDim,
1688 const t_state& state,
1689 const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
1690 const gmx_ekindata_t& ekind)
1692 GMX_ASSERT(inputRecord.eI == IntegrationAlgorithm::MD || inputRecord.eI == IntegrationAlgorithm::SD1,
1693 "Only leap-frog is supported here");
1695 // Cast to real for faster code, no loss in precision
1696 const real dt = inputRecord.delta_t;
1698 const int nth = gmx_omp_nthreads_get(ModuleMultiThread::Update);
1700 #pragma omp parallel for num_threads(nth) schedule(static)
1701 for (int th = 0; th < nth; th++)
1705 int start_th, end_th;
1706 getThreadAtomRange(nth, th, homenr, &start_th, &end_th);
1708 const rvec* x_rvec = state.x.rvec_array();
1709 rvec* xp_rvec = xp_.rvec_array();
1710 const rvec* v_rvec = state.v.rvec_array();
1711 const rvec* f_rvec = as_rvec_array(f.unpaddedConstArrayRef().data());
1713 doUpdateMDDoNotUpdateVelocities(
1714 start_th, end_th, dt, x_rvec, xp_rvec, v_rvec, f_rvec, havePartiallyFrozenAtoms, invmass, invMassPerDim, ekind);
1716 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR