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 the velocities of virtual sites to zero */
259 static void clearVsiteVelocities(int start, int nrend, const unsigned short* particleType, rvec* gmx_restrict v)
261 for (int a = start; a < nrend; a++)
263 if (particleType[a] == eptVSite)
270 /*! \brief Sets whether we store the updated velocities */
271 enum class StoreUpdatedVelocities
273 yes, //!< Store the updated velocities
274 no //!< Do not store the updated velocities
277 /*! \brief Sets the number of different temperature coupling values */
278 enum class NumTempScaleValues
280 single, //!< Single T-scaling value (either one group or all values =1)
281 multiple //!< Multiple T-scaling values, need to use T-group indices
284 /*! \brief Sets if to apply no or diagonal Parrinello-Rahman pressure scaling
286 * Note that this enum is only used in updateMDLeapfrogSimple(), which does
287 * not handle fully anistropic Parrinello-Rahman scaling, so we only have
288 * options \p no and \p diagonal here and no anistropic option.
290 enum class ApplyParrinelloRahmanVScaling
292 no, //!< Do not apply velocity scaling (not a PR-coupling run or step)
293 diagonal //!< Apply velocity scaling using a diagonal matrix
296 /*! \brief Integrate using leap-frog with T-scaling and optionally diagonal Parrinello-Rahman p-coupling
298 * \tparam storeUpdatedVelocities Tells whether we should store the updated velocities
299 * \tparam numTempScaleValues The number of different T-couple values
300 * \tparam applyPRVScaling Apply Parrinello-Rahman velocity scaling
301 * \param[in] start Index of first atom to update
302 * \param[in] nrend Last atom to update: \p nrend - 1
303 * \param[in] dt The time step
304 * \param[in] dtPressureCouple Time step for pressure coupling
305 * \param[in] invMassPerDim 1/mass per atom and dimension
306 * \param[in] tcstat Temperature coupling information
307 * \param[in] cTC T-coupling group index per atom
308 * \param[in] pRVScaleMatrixDiagonal Parrinello-Rahman v-scale matrix diagonal
309 * \param[in] x Input coordinates
310 * \param[out] xprime Updated coordinates
311 * \param[inout] v Velocities
312 * \param[in] f Forces
314 * We expect this template to get good SIMD acceleration by most compilers,
315 * unlike the more complex general template.
316 * Note that we might get even better SIMD acceleration when we introduce
317 * aligned (and padded) memory, possibly with some hints for the compilers.
319 template<StoreUpdatedVelocities storeUpdatedVelocities, NumTempScaleValues numTempScaleValues, ApplyParrinelloRahmanVScaling applyPRVScaling>
320 static void updateMDLeapfrogSimple(int start,
323 real dtPressureCouple,
324 const rvec* gmx_restrict invMassPerDim,
325 gmx::ArrayRef<const t_grp_tcstat> tcstat,
326 const unsigned short* cTC,
327 const rvec pRVScaleMatrixDiagonal,
328 const rvec* gmx_restrict x,
329 rvec* gmx_restrict xprime,
330 rvec* gmx_restrict v,
331 const rvec* gmx_restrict f)
335 if (numTempScaleValues == NumTempScaleValues::single)
337 lambdaGroup = tcstat[0].lambda;
340 for (int a = start; a < nrend; a++)
342 if (numTempScaleValues == NumTempScaleValues::multiple)
344 lambdaGroup = tcstat[cTC[a]].lambda;
347 for (int d = 0; d < DIM; d++)
349 /* Note that using rvec invMassPerDim results in more efficient
350 * SIMD code, but this increases the cache pressure.
351 * For large systems with PME on the CPU this slows down the
352 * (then already slow) update by 20%. If all data remains in cache,
353 * using rvec is much faster.
355 real vNew = lambdaGroup * v[a][d] + f[a][d] * invMassPerDim[a][d] * dt;
357 if (applyPRVScaling == ApplyParrinelloRahmanVScaling::diagonal)
359 vNew -= dtPressureCouple * pRVScaleMatrixDiagonal[d] * v[a][d];
361 if (storeUpdatedVelocities == StoreUpdatedVelocities::yes)
365 xprime[a][d] = x[a][d] + vNew * dt;
370 #if GMX_SIMD && GMX_SIMD_HAVE_REAL
371 # define GMX_HAVE_SIMD_UPDATE 1
373 # define GMX_HAVE_SIMD_UPDATE 0
376 #if GMX_HAVE_SIMD_UPDATE
378 /*! \brief Load (aligned) the contents of GMX_SIMD_REAL_WIDTH rvec elements sequentially into 3 SIMD registers
380 * The loaded output is:
381 * \p r0: { r[index][XX], r[index][YY], ... }
383 * \p r2: { ..., r[index+GMX_SIMD_REAL_WIDTH-1][YY], r[index+GMX_SIMD_REAL_WIDTH-1][ZZ] }
385 * \param[in] r Real to an rvec array, has to be aligned to SIMD register width
386 * \param[in] index Index of the first rvec triplet of reals to load
387 * \param[out] r0 Pointer to first SIMD register
388 * \param[out] r1 Pointer to second SIMD register
389 * \param[out] r2 Pointer to third SIMD register
391 static inline void simdLoadRvecs(const rvec* r, int index, SimdReal* r0, SimdReal* r1, SimdReal* r2)
393 const real* realPtr = r[index];
395 GMX_ASSERT(isSimdAligned(realPtr), "Pointer should be SIMD aligned");
397 *r0 = simdLoad(realPtr + 0 * GMX_SIMD_REAL_WIDTH);
398 *r1 = simdLoad(realPtr + 1 * GMX_SIMD_REAL_WIDTH);
399 *r2 = simdLoad(realPtr + 2 * GMX_SIMD_REAL_WIDTH);
402 /*! \brief Store (aligned) 3 SIMD registers sequentially to GMX_SIMD_REAL_WIDTH rvec elements
404 * The stored output is:
405 * \p r[index] = { { r0[0], r0[1], ... }
407 * \p r[index+GMX_SIMD_REAL_WIDTH-1] = { ... , r2[GMX_SIMD_REAL_WIDTH-2], r2[GMX_SIMD_REAL_WIDTH-1] }
409 * \param[out] r Pointer to an rvec array, has to be aligned to SIMD register width
410 * \param[in] index Index of the first rvec triplet of reals to store to
411 * \param[in] r0 First SIMD register
412 * \param[in] r1 Second SIMD register
413 * \param[in] r2 Third SIMD register
415 static inline void simdStoreRvecs(rvec* r, int index, SimdReal r0, SimdReal r1, SimdReal r2)
417 real* realPtr = r[index];
419 GMX_ASSERT(isSimdAligned(realPtr), "Pointer should be SIMD aligned");
421 store(realPtr + 0 * GMX_SIMD_REAL_WIDTH, r0);
422 store(realPtr + 1 * GMX_SIMD_REAL_WIDTH, r1);
423 store(realPtr + 2 * GMX_SIMD_REAL_WIDTH, r2);
426 /*! \brief Integrate using leap-frog with single group T-scaling and SIMD
428 * \tparam storeUpdatedVelocities Tells whether we should store the updated velocities
429 * \param[in] start Index of first atom to update
430 * \param[in] nrend Last atom to update: \p nrend - 1
431 * \param[in] dt The time step
432 * \param[in] invMass 1/mass per atom
433 * \param[in] tcstat Temperature coupling information
434 * \param[in] x Input coordinates
435 * \param[out] xprime Updated coordinates
436 * \param[inout] v Velocities
437 * \param[in] f Forces
439 template<StoreUpdatedVelocities storeUpdatedVelocities>
440 static void updateMDLeapfrogSimpleSimd(int start,
443 const real* gmx_restrict invMass,
444 gmx::ArrayRef<const t_grp_tcstat> tcstat,
445 const rvec* gmx_restrict x,
446 rvec* gmx_restrict xprime,
447 rvec* gmx_restrict v,
448 const rvec* gmx_restrict f)
450 SimdReal timestep(dt);
451 SimdReal lambdaSystem(tcstat[0].lambda);
453 /* We declare variables here, since code is often slower when declaring them inside the loop */
455 /* Note: We should implement a proper PaddedVector, so we don't need this check */
456 GMX_ASSERT(isSimdAligned(invMass), "invMass should be aligned");
458 for (int a = start; a < nrend; a += GMX_SIMD_REAL_WIDTH)
460 SimdReal invMass0, invMass1, invMass2;
461 expandScalarsToTriplets(simdLoad(invMass + a), &invMass0, &invMass1, &invMass2);
465 simdLoadRvecs(v, a, &v0, &v1, &v2);
466 simdLoadRvecs(f, a, &f0, &f1, &f2);
468 v0 = fma(f0 * invMass0, timestep, lambdaSystem * v0);
469 v1 = fma(f1 * invMass1, timestep, lambdaSystem * v1);
470 v2 = fma(f2 * invMass2, timestep, lambdaSystem * v2);
472 if (storeUpdatedVelocities == StoreUpdatedVelocities::yes)
474 simdStoreRvecs(v, a, v0, v1, v2);
478 simdLoadRvecs(x, a, &x0, &x1, &x2);
480 SimdReal xprime0 = fma(v0, timestep, x0);
481 SimdReal xprime1 = fma(v1, timestep, x1);
482 SimdReal xprime2 = fma(v2, timestep, x2);
484 simdStoreRvecs(xprime, a, xprime0, xprime1, xprime2);
488 #endif // GMX_HAVE_SIMD_UPDATE
490 /*! \brief Sets the NEMD acceleration type */
491 enum class AccelerationType
497 /*! \brief Integrate using leap-frog with support for everything.
499 * \tparam accelerationType Type of NEMD acceleration.
500 * \param[in] start Index of first atom to update.
501 * \param[in] nrend Last atom to update: \p nrend - 1.
502 * \param[in] doNoseHoover If to apply Nose-Hoover T-coupling.
503 * \param[in] dt The time step.
504 * \param[in] dtPressureCouple Time step for pressure coupling, is 0 when no pressure
505 * coupling should be applied at this step.
506 * \param[in] md Atom properties.
507 * \param[in] ekind Kinetic energy data.
508 * \param[in] box The box dimensions.
509 * \param[in] x Input coordinates.
510 * \param[out] xprime Updated coordinates.
511 * \param[inout] v Velocities.
512 * \param[in] f Forces.
513 * \param[in] nh_vxi Nose-Hoover velocity scaling factors.
514 * \param[in] nsttcouple Frequency of the temperature coupling steps.
515 * \param[in] M Parrinello-Rahman scaling matrix.
517 template<AccelerationType accelerationType>
518 static void updateMDLeapfrogGeneral(int start,
522 real dtPressureCouple,
524 const gmx_ekindata_t* ekind,
526 const rvec* gmx_restrict x,
527 rvec* gmx_restrict xprime,
528 rvec* gmx_restrict v,
529 const rvec* gmx_restrict f,
530 const double* gmx_restrict nh_vxi,
531 const int nsttcouple,
534 /* This is a version of the leap-frog integrator that supports
535 * all combinations of T-coupling, P-coupling and NEMD.
536 * Nose-Hoover uses the reversible leap-frog integrator from
537 * Holian et al. Phys Rev E 52(3) : 2338, 1995
540 gmx::ArrayRef<const t_grp_tcstat> tcstat = ekind->tcstat;
541 const unsigned short* cTC = md->cTC;
543 const rvec* gmx_restrict invMassPerDim = md->invMassPerDim;
545 /* Initialize group values, changed later when multiple groups are used */
548 real omega_Z = 2 * static_cast<real>(M_PI) / box[ZZ][ZZ];
550 for (int n = start; n < nrend; n++)
556 real lg = tcstat[gt].lambda;
559 real cosineZ, vCosine;
560 #ifdef __INTEL_COMPILER
561 # pragma warning(disable : 280)
563 switch (accelerationType)
565 case AccelerationType::none: copy_rvec(v[n], vRel); break;
566 case AccelerationType::cosine:
567 cosineZ = std::cos(x[n][ZZ] * omega_Z);
568 vCosine = cosineZ * ekind->cosacc.vcos;
569 /* Avoid scaling the cosine profile velocity */
570 copy_rvec(v[n], vRel);
578 /* Here we account for multiple time stepping, by increasing
579 * the Nose-Hoover correction by nsttcouple
580 * TODO: This can be pre-computed.
582 factorNH = 0.5 * nsttcouple * dt * nh_vxi[gt];
585 for (int d = 0; d < DIM; d++)
587 real vNew = (lg * vRel[d]
588 + (f[n][d] * invMassPerDim[n][d] * dt - factorNH * vRel[d]
589 - dtPressureCouple * iprod(M[d], vRel)))
591 switch (accelerationType)
593 case AccelerationType::none: break;
594 case AccelerationType::cosine:
597 /* Add back the mean velocity and apply acceleration */
598 vNew += vCosine + cosineZ * ekind->cosacc.cos_accel * dt;
603 xprime[n][d] = x[n][d] + vNew * dt;
608 /*! \brief Handles the Leap-frog MD x and v integration */
609 static void do_update_md(int start,
613 const rvec* gmx_restrict x,
614 rvec* gmx_restrict xprime,
615 rvec* gmx_restrict v,
616 const rvec* gmx_restrict f,
619 const int nsttcouple,
620 const int nstpcouple,
622 const gmx_ekindata_t* ekind,
624 const double* gmx_restrict nh_vxi,
627 GMX_ASSERT(nrend == start || xprime != x,
628 "For SIMD optimization certain compilers need to have xprime != x");
630 /* Note: Berendsen pressure scaling is handled after do_update_md() */
631 bool doTempCouple = (etc != etcNO && do_per_step(step + nsttcouple - 1, nsttcouple));
632 bool doNoseHoover = (etc == etcNOSEHOOVER && doTempCouple);
633 bool doParrinelloRahman =
634 (epc == epcPARRINELLORAHMAN && do_per_step(step + nstpcouple - 1, nstpcouple));
635 bool doPROffDiagonal = (doParrinelloRahman && (M[YY][XX] != 0 || M[ZZ][XX] != 0 || M[ZZ][YY] != 0));
637 real dtPressureCouple = (doParrinelloRahman ? nstpcouple * dt : 0);
639 /* NEMD cosine acceleration is applied in updateMDLeapFrogGeneral */
640 bool doAcceleration = (ekind->cosacc.cos_accel != 0);
642 if (doNoseHoover || doPROffDiagonal || doAcceleration)
645 if (!doParrinelloRahman)
647 /* We should not apply PR scaling at this step */
657 updateMDLeapfrogGeneral<AccelerationType::none>(
658 start, nrend, doNoseHoover, dt, dtPressureCouple, md, ekind, box, x, xprime, v, f, nh_vxi, nsttcouple, stepM);
662 updateMDLeapfrogGeneral<AccelerationType::cosine>(
663 start, nrend, doNoseHoover, dt, dtPressureCouple, md, ekind, box, x, xprime, v, f, nh_vxi, nsttcouple, stepM);
668 /* Use a simple and thus more efficient integration loop. */
669 /* The simple loop does not check for particle type (so it can
670 * be vectorized), which means we need to clear the velocities
671 * of virtual sites in advance, when present. Note that vsite
672 * velocities are computed after update and constraints from
673 * their displacement.
677 /* Note: The overhead of this loop is completely neligible */
678 clearVsiteVelocities(start, nrend, md->ptype, v);
681 /* We determine if we have a single T-coupling lambda value for all
682 * atoms. That allows for better SIMD acceleration in the template.
683 * If we do not do temperature coupling (in the run or this step),
684 * all scaling values are 1, so we effectively have a single value.
686 bool haveSingleTempScaleValue = (!doTempCouple || ekind->ngtc == 1);
688 /* Extract some pointers needed by all cases */
689 const unsigned short* cTC = md->cTC;
690 gmx::ArrayRef<const t_grp_tcstat> tcstat = ekind->tcstat;
691 const rvec* invMassPerDim = md->invMassPerDim;
693 if (doParrinelloRahman)
695 GMX_ASSERT(!doPROffDiagonal,
696 "updateMDLeapfrogSimple only support diagonal Parrinello-Rahman scaling "
700 for (int d = 0; d < DIM; d++)
705 if (haveSingleTempScaleValue)
707 updateMDLeapfrogSimple<StoreUpdatedVelocities::yes, NumTempScaleValues::single, ApplyParrinelloRahmanVScaling::diagonal>(
708 start, nrend, dt, dtPressureCouple, invMassPerDim, tcstat, cTC, diagM, x, xprime, v, f);
712 updateMDLeapfrogSimple<StoreUpdatedVelocities::yes, NumTempScaleValues::multiple, ApplyParrinelloRahmanVScaling::diagonal>(
713 start, nrend, dt, dtPressureCouple, invMassPerDim, tcstat, cTC, diagM, x, xprime, v, f);
718 if (haveSingleTempScaleValue)
720 /* Note that modern compilers are pretty good at vectorizing
721 * updateMDLeapfrogSimple(). But the SIMD version will still
722 * be faster because invMass lowers the cache pressure
723 * compared to invMassPerDim.
725 #if GMX_HAVE_SIMD_UPDATE
726 /* Check if we can use invmass instead of invMassPerDim */
727 if (!md->havePartiallyFrozenAtoms)
729 updateMDLeapfrogSimpleSimd<StoreUpdatedVelocities::yes>(
730 start, nrend, dt, md->invmass, tcstat, x, xprime, v, f);
735 updateMDLeapfrogSimple<StoreUpdatedVelocities::yes, NumTempScaleValues::single, ApplyParrinelloRahmanVScaling::no>(
736 start, nrend, dt, dtPressureCouple, invMassPerDim, tcstat, cTC, nullptr, x, xprime, v, f);
741 updateMDLeapfrogSimple<StoreUpdatedVelocities::yes, NumTempScaleValues::multiple, ApplyParrinelloRahmanVScaling::no>(
742 start, nrend, dt, dtPressureCouple, invMassPerDim, tcstat, cTC, nullptr, x, xprime, v, f);
747 /*! \brief Handles the Leap-frog MD x and v integration */
748 static void doUpdateMDDoNotUpdateVelocities(int start,
751 const rvec* gmx_restrict x,
752 rvec* gmx_restrict xprime,
753 rvec* gmx_restrict v,
754 const rvec* gmx_restrict f,
756 const gmx_ekindata_t& ekind)
758 GMX_ASSERT(nrend == start || xprime != x,
759 "For SIMD optimization certain compilers need to have xprime != x");
761 gmx::ArrayRef<const t_grp_tcstat> tcstat = ekind.tcstat;
763 /* Check if we can use invmass instead of invMassPerDim */
764 #if GMX_HAVE_SIMD_UPDATE
765 if (!md.havePartiallyFrozenAtoms)
767 updateMDLeapfrogSimpleSimd<StoreUpdatedVelocities::no>(
768 start, nrend, dt, md.invmass, tcstat, x, xprime, v, f);
773 updateMDLeapfrogSimple<StoreUpdatedVelocities::no, NumTempScaleValues::single, ApplyParrinelloRahmanVScaling::no>(
774 start, nrend, dt, dt, md.invMassPerDim, tcstat, nullptr, nullptr, x, xprime, v, f);
778 static void do_update_vv_vel(int start,
781 const ivec nFreeze[],
782 const real invmass[],
783 const unsigned short ptype[],
784 const unsigned short cFREEZE[],
797 g = 0.25 * dt * veta * alpha;
799 mv2 = gmx::series_sinhx(g);
806 for (n = start; n < nrend; n++)
808 real w_dt = invmass[n] * dt;
814 for (d = 0; d < DIM; d++)
816 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
818 v[n][d] = mv1 * (mv1 * v[n][d] + 0.5 * (w_dt * mv2 * f[n][d]));
826 } /* do_update_vv_vel */
828 static void do_update_vv_pos(int start,
831 const ivec nFreeze[],
832 const unsigned short ptype[],
833 const unsigned short cFREEZE[],
844 /* Would it make more sense if Parrinello-Rahman was put here? */
849 mr2 = gmx::series_sinhx(g);
857 for (n = start; n < nrend; n++)
865 for (d = 0; d < DIM; d++)
867 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
869 xprime[n][d] = mr1 * (mr1 * x[n][d] + mr2 * dt * v[n][d]);
873 xprime[n][d] = x[n][d];
877 } /* do_update_vv_pos */
879 gmx_stochd_t::gmx_stochd_t(const t_inputrec& inputRecord)
881 const t_grpopts* opts = &inputRecord.opts;
882 const int ngtc = opts->ngtc;
884 if (inputRecord.eI == eiBD)
888 else if (EI_SD(inputRecord.eI))
893 for (int gt = 0; gt < ngtc; gt++)
895 if (opts->tau_t[gt] > 0)
897 sdc[gt].em = std::exp(-inputRecord.delta_t / opts->tau_t[gt]);
901 /* No friction and noise on this group */
906 else if (ETC_ANDERSEN(inputRecord.etc))
908 randomize_group.resize(ngtc);
909 boltzfac.resize(ngtc);
911 /* for now, assume that all groups, if randomized, are randomized at the same rate, i.e. tau_t is the same. */
912 /* since constraint groups don't necessarily match up with temperature groups! This is checked in readir.c */
914 for (int gt = 0; gt < ngtc; gt++)
916 real reft = std::max<real>(0, opts->ref_t[gt]);
917 if ((opts->tau_t[gt] > 0)
918 && (reft > 0)) /* tau_t or ref_t = 0 means that no randomization is done */
920 randomize_group[gt] = true;
921 boltzfac[gt] = BOLTZ * opts->ref_t[gt];
925 randomize_group[gt] = false;
931 void Update::Impl::update_temperature_constants(const t_inputrec& inputRecord)
933 if (inputRecord.eI == eiBD)
935 if (inputRecord.bd_fric != 0)
937 for (int gt = 0; gt < inputRecord.opts.ngtc; gt++)
939 sd_.bd_rf[gt] = std::sqrt(2.0 * BOLTZ * inputRecord.opts.ref_t[gt]
940 / (inputRecord.bd_fric * inputRecord.delta_t));
945 for (int gt = 0; gt < inputRecord.opts.ngtc; gt++)
947 sd_.bd_rf[gt] = std::sqrt(2.0 * BOLTZ * inputRecord.opts.ref_t[gt]);
951 if (inputRecord.eI == eiSD1)
953 for (int gt = 0; gt < inputRecord.opts.ngtc; gt++)
955 real kT = BOLTZ * inputRecord.opts.ref_t[gt];
956 /* The mass is accounted for later, since this differs per atom */
957 sd_.sdsig[gt].V = std::sqrt(kT * (1 - sd_.sdc[gt].em * sd_.sdc[gt].em));
962 Update::Impl::Impl(const t_inputrec& inputRecord, BoxDeformation* boxDeformation) :
964 deform_(boxDeformation)
966 update_temperature_constants(inputRecord);
967 xp_.resizeWithPadding(0);
970 void Update::setNumAtoms(int numAtoms)
973 impl_->xp()->resizeWithPadding(numAtoms);
976 /*! \brief Sets the SD update type */
977 enum class SDUpdate : int
980 FrictionAndNoiseOnly,
984 /*! \brief SD integrator update
986 * Two phases are required in the general case of a constrained
987 * update, the first phase from the contribution of forces, before
988 * applying constraints, and then a second phase applying the friction
989 * and noise, and then further constraining. For details, see
992 * Without constraints, the two phases can be combined, for
995 * Thus three instantiations of this templated function will be made,
996 * two with only one contribution, and one with both contributions. */
997 template<SDUpdate updateType>
998 static void doSDUpdateGeneral(const gmx_stochd_t& sd,
1002 const ivec nFreeze[],
1003 const real invmass[],
1004 const unsigned short ptype[],
1005 const unsigned short cFREEZE[],
1006 const unsigned short cTC[],
1013 const int* gatindex)
1015 // cTC and cFREEZE can be nullptr any time, but various
1016 // instantiations do not make sense with particular pointer
1018 if (updateType == SDUpdate::ForcesOnly)
1020 GMX_ASSERT(f != nullptr, "SD update with only forces requires forces");
1021 GMX_ASSERT(cTC == nullptr, "SD update with only forces cannot handle temperature groups");
1023 if (updateType == SDUpdate::FrictionAndNoiseOnly)
1025 GMX_ASSERT(f == nullptr, "SD update with only noise cannot handle forces");
1027 if (updateType == SDUpdate::Combined)
1029 GMX_ASSERT(f != nullptr, "SD update with forces and noise requires forces");
1032 // Even 0 bits internal counter gives 2x64 ints (more than enough for three table lookups)
1033 gmx::ThreeFry2x64<0> rng(seed, gmx::RandomDomain::UpdateCoordinates);
1034 gmx::TabulatedNormalDistribution<real, 14> dist;
1036 for (int n = start; n < nrend; n++)
1038 int globalAtomIndex = gatindex ? gatindex[n] : n;
1039 rng.restart(step, globalAtomIndex);
1042 real inverseMass = invmass[n];
1043 real invsqrtMass = std::sqrt(inverseMass);
1045 int freezeGroup = cFREEZE ? cFREEZE[n] : 0;
1046 int temperatureGroup = cTC ? cTC[n] : 0;
1048 for (int d = 0; d < DIM; d++)
1050 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[freezeGroup][d])
1052 if (updateType == SDUpdate::ForcesOnly)
1054 real vn = v[n][d] + inverseMass * f[n][d] * dt;
1056 // Simple position update.
1057 xprime[n][d] = x[n][d] + v[n][d] * dt;
1059 else if (updateType == SDUpdate::FrictionAndNoiseOnly)
1062 v[n][d] = (vn * sd.sdc[temperatureGroup].em
1063 + invsqrtMass * sd.sdsig[temperatureGroup].V * dist(rng));
1064 // The previous phase already updated the
1065 // positions with a full v*dt term that must
1066 // now be half removed.
1067 xprime[n][d] = xprime[n][d] + 0.5 * (v[n][d] - vn) * dt;
1071 real vn = v[n][d] + inverseMass * f[n][d] * dt;
1072 v[n][d] = (vn * sd.sdc[temperatureGroup].em
1073 + invsqrtMass * sd.sdsig[temperatureGroup].V * dist(rng));
1074 // Here we include half of the friction+noise
1075 // update of v into the position update.
1076 xprime[n][d] = x[n][d] + 0.5 * (vn + v[n][d]) * dt;
1081 // When using constraints, the update is split into
1082 // two phases, but we only need to zero the update of
1083 // virtual, shell or frozen particles in at most one
1085 if (updateType != SDUpdate::FrictionAndNoiseOnly)
1088 xprime[n][d] = x[n][d];
1095 static void do_update_sd(int start,
1099 const rvec* gmx_restrict x,
1100 rvec* gmx_restrict xprime,
1101 rvec* gmx_restrict v,
1102 const rvec* gmx_restrict f,
1103 const ivec nFreeze[],
1104 const real invmass[],
1105 const unsigned short ptype[],
1106 const unsigned short cFREEZE[],
1107 const unsigned short cTC[],
1109 const t_commrec* cr,
1110 const gmx_stochd_t& sd,
1111 bool haveConstraints)
1113 if (haveConstraints)
1115 // With constraints, the SD update is done in 2 parts
1116 doSDUpdateGeneral<SDUpdate::ForcesOnly>(
1117 sd, start, nrend, dt, nFreeze, invmass, ptype, cFREEZE, nullptr, x, xprime, v, f, step, seed, nullptr);
1121 doSDUpdateGeneral<SDUpdate::Combined>(sd,
1136 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1140 static void do_update_bd(int start,
1144 const rvec* gmx_restrict x,
1145 rvec* gmx_restrict xprime,
1146 rvec* gmx_restrict v,
1147 const rvec* gmx_restrict f,
1148 const ivec nFreeze[],
1149 const real invmass[],
1150 const unsigned short ptype[],
1151 const unsigned short cFREEZE[],
1152 const unsigned short cTC[],
1153 real friction_coefficient,
1156 const int* gatindex)
1158 /* note -- these appear to be full step velocities . . . */
1163 // Use 1 bit of internal counters to give us 2*2 64-bits values per stream
1164 // Each 64-bit value is enough for 4 normal distribution table numbers.
1165 gmx::ThreeFry2x64<0> rng(seed, gmx::RandomDomain::UpdateCoordinates);
1166 gmx::TabulatedNormalDistribution<real, 14> dist;
1168 if (friction_coefficient != 0)
1170 invfr = 1.0 / friction_coefficient;
1173 for (n = start; (n < nrend); n++)
1175 int ng = gatindex ? gatindex[n] : n;
1177 rng.restart(step, ng);
1188 for (d = 0; (d < DIM); d++)
1190 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
1192 if (friction_coefficient != 0)
1194 vn = invfr * f[n][d] + rf[gt] * dist(rng);
1198 /* NOTE: invmass = 2/(mass*friction_constant*dt) */
1199 vn = 0.5 * invmass[n] * f[n][d] * dt
1200 + std::sqrt(0.5 * invmass[n]) * rf[gt] * dist(rng);
1204 xprime[n][d] = x[n][d] + vn * dt;
1209 xprime[n][d] = x[n][d];
1215 extern void init_ekinstate(ekinstate_t* ekinstate, const t_inputrec* ir)
1217 ekinstate->ekin_n = ir->opts.ngtc;
1218 snew(ekinstate->ekinh, ekinstate->ekin_n);
1219 snew(ekinstate->ekinf, ekinstate->ekin_n);
1220 snew(ekinstate->ekinh_old, ekinstate->ekin_n);
1221 ekinstate->ekinscalef_nhc.resize(ekinstate->ekin_n);
1222 ekinstate->ekinscaleh_nhc.resize(ekinstate->ekin_n);
1223 ekinstate->vscale_nhc.resize(ekinstate->ekin_n);
1224 ekinstate->dekindl = 0;
1225 ekinstate->mvcos = 0;
1226 ekinstate->hasReadEkinState = false;
1229 void update_ekinstate(ekinstate_t* ekinstate, const gmx_ekindata_t* ekind)
1233 for (i = 0; i < ekinstate->ekin_n; i++)
1235 copy_mat(ekind->tcstat[i].ekinh, ekinstate->ekinh[i]);
1236 copy_mat(ekind->tcstat[i].ekinf, ekinstate->ekinf[i]);
1237 copy_mat(ekind->tcstat[i].ekinh_old, ekinstate->ekinh_old[i]);
1238 ekinstate->ekinscalef_nhc[i] = ekind->tcstat[i].ekinscalef_nhc;
1239 ekinstate->ekinscaleh_nhc[i] = ekind->tcstat[i].ekinscaleh_nhc;
1240 ekinstate->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
1243 copy_mat(ekind->ekin, ekinstate->ekin_total);
1244 ekinstate->dekindl = ekind->dekindl;
1245 ekinstate->mvcos = ekind->cosacc.mvcos;
1248 void restore_ekinstate_from_state(const t_commrec* cr, gmx_ekindata_t* ekind, const ekinstate_t* ekinstate)
1254 for (i = 0; i < ekinstate->ekin_n; i++)
1256 copy_mat(ekinstate->ekinh[i], ekind->tcstat[i].ekinh);
1257 copy_mat(ekinstate->ekinf[i], ekind->tcstat[i].ekinf);
1258 copy_mat(ekinstate->ekinh_old[i], ekind->tcstat[i].ekinh_old);
1259 ekind->tcstat[i].ekinscalef_nhc = ekinstate->ekinscalef_nhc[i];
1260 ekind->tcstat[i].ekinscaleh_nhc = ekinstate->ekinscaleh_nhc[i];
1261 ekind->tcstat[i].vscale_nhc = ekinstate->vscale_nhc[i];
1264 copy_mat(ekinstate->ekin_total, ekind->ekin);
1266 ekind->dekindl = ekinstate->dekindl;
1267 ekind->cosacc.mvcos = ekinstate->mvcos;
1268 n = ekinstate->ekin_n;
1273 gmx_bcast(sizeof(n), &n, cr->mpi_comm_mygroup);
1274 for (i = 0; i < n; i++)
1276 gmx_bcast(DIM * DIM * sizeof(ekind->tcstat[i].ekinh[0][0]),
1277 ekind->tcstat[i].ekinh[0],
1278 cr->mpi_comm_mygroup);
1279 gmx_bcast(DIM * DIM * sizeof(ekind->tcstat[i].ekinf[0][0]),
1280 ekind->tcstat[i].ekinf[0],
1281 cr->mpi_comm_mygroup);
1282 gmx_bcast(DIM * DIM * sizeof(ekind->tcstat[i].ekinh_old[0][0]),
1283 ekind->tcstat[i].ekinh_old[0],
1284 cr->mpi_comm_mygroup);
1286 gmx_bcast(sizeof(ekind->tcstat[i].ekinscalef_nhc),
1287 &(ekind->tcstat[i].ekinscalef_nhc),
1288 cr->mpi_comm_mygroup);
1289 gmx_bcast(sizeof(ekind->tcstat[i].ekinscaleh_nhc),
1290 &(ekind->tcstat[i].ekinscaleh_nhc),
1291 cr->mpi_comm_mygroup);
1292 gmx_bcast(sizeof(ekind->tcstat[i].vscale_nhc), &(ekind->tcstat[i].vscale_nhc), cr->mpi_comm_mygroup);
1294 gmx_bcast(DIM * DIM * sizeof(ekind->ekin[0][0]), ekind->ekin[0], cr->mpi_comm_mygroup);
1296 gmx_bcast(sizeof(ekind->dekindl), &ekind->dekindl, cr->mpi_comm_mygroup);
1297 gmx_bcast(sizeof(ekind->cosacc.mvcos), &ekind->cosacc.mvcos, cr->mpi_comm_mygroup);
1301 void getThreadAtomRange(int numThreads, int threadIndex, int numAtoms, int* startAtom, int* endAtom)
1303 #if GMX_HAVE_SIMD_UPDATE
1304 constexpr int blockSize = GMX_SIMD_REAL_WIDTH;
1306 constexpr int blockSize = 1;
1308 int numBlocks = (numAtoms + blockSize - 1) / blockSize;
1310 *startAtom = ((numBlocks * threadIndex) / numThreads) * blockSize;
1311 *endAtom = ((numBlocks * (threadIndex + 1)) / numThreads) * blockSize;
1312 if (threadIndex == numThreads - 1)
1314 *endAtom = numAtoms;
1318 void Update::Impl::update_sd_second_half(const t_inputrec& inputRecord,
1321 const t_mdatoms* md,
1323 const t_commrec* cr,
1325 gmx_wallcycle_t wcycle,
1326 gmx::Constraints* constr,
1334 if (inputRecord.eI == eiSD1)
1336 int homenr = md->homenr;
1338 /* Cast delta_t from double to real to make the integrators faster.
1339 * The only reason for having delta_t double is to get accurate values
1340 * for t=delta_t*step when step is larger than float precision.
1341 * For integration dt the accuracy of real suffices, since with
1342 * integral += dt*integrand the increment is nearly always (much) smaller
1343 * than the integral (and the integrand has real precision).
1345 real dt = inputRecord.delta_t;
1347 wallcycle_start(wcycle, ewcUPDATE);
1349 int nth = gmx_omp_nthreads_get(emntUpdate);
1351 #pragma omp parallel for num_threads(nth) schedule(static)
1352 for (int th = 0; th < nth; th++)
1356 int start_th, end_th;
1357 getThreadAtomRange(nth, th, homenr, &start_th, &end_th);
1359 doSDUpdateGeneral<SDUpdate::FrictionAndNoiseOnly>(
1364 inputRecord.opts.nFreeze,
1369 state->x.rvec_array(),
1371 state->v.rvec_array(),
1374 inputRecord.ld_seed,
1375 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1377 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1379 inc_nrnb(nrnb, eNR_UPDATE, homenr);
1380 wallcycle_stop(wcycle, ewcUPDATE);
1382 /* Constrain the coordinates upd->xp for half a time step */
1383 bool computeVirial = false;
1384 constr->apply(do_log,
1389 state->x.arrayRefWithPadding(),
1390 xp_.arrayRefWithPadding(),
1393 state->lambda[efptBONDED],
1395 state->v.arrayRefWithPadding(),
1398 ConstraintVariable::Positions);
1402 void Update::Impl::finish_update(const t_inputrec& inputRecord,
1403 const t_mdatoms* md,
1405 gmx_wallcycle_t wcycle,
1406 const bool haveConstraints)
1408 /* NOTE: Currently we always integrate to a temporary buffer and
1409 * then copy the results back here.
1412 wallcycle_start_nocount(wcycle, ewcUPDATE);
1414 const int homenr = md->homenr;
1415 auto xp = makeConstArrayRef(xp_).subArray(0, homenr);
1416 auto x = makeArrayRef(state->x).subArray(0, homenr);
1418 if (md->havePartiallyFrozenAtoms && haveConstraints)
1420 /* We have atoms that are frozen along some, but not all dimensions,
1421 * then constraints will have moved them also along the frozen dimensions.
1422 * To freeze such degrees of freedom we do not copy them back here.
1424 const ivec* nFreeze = inputRecord.opts.nFreeze;
1426 for (int i = 0; i < homenr; i++)
1428 const int g = md->cFREEZE[i];
1430 for (int d = 0; d < DIM; d++)
1432 if (nFreeze[g][d] == 0)
1441 /* We have no frozen atoms or fully frozen atoms which have not
1442 * been moved by the update, so we can simply copy all coordinates.
1444 int gmx_unused nth = gmx_omp_nthreads_get(emntUpdate);
1445 #pragma omp parallel for num_threads(nth) schedule(static)
1446 for (int i = 0; i < homenr; i++)
1448 // Trivial statement, does not throw
1453 wallcycle_stop(wcycle, ewcUPDATE);
1456 void Update::Impl::update_coords(const t_inputrec& inputRecord,
1458 const t_mdatoms* md,
1460 const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
1461 const t_fcdata& fcdata,
1462 const gmx_ekindata_t* ekind,
1465 const t_commrec* cr,
1466 const bool haveConstraints)
1468 /* Running the velocity half does nothing except for velocity verlet */
1469 if ((updatePart == etrtVELOCITY1 || updatePart == etrtVELOCITY2) && !EI_VV(inputRecord.eI))
1471 gmx_incons("update_coords called for velocity without VV integrator");
1474 int homenr = md->homenr;
1476 /* Cast to real for faster code, no loss in precision (see comment above) */
1477 real dt = inputRecord.delta_t;
1479 /* We need to update the NMR restraint history when time averaging is used */
1480 if (state->flags & (1 << estDISRE_RM3TAV))
1482 update_disres_history(*fcdata.disres, &state->hist);
1484 if (state->flags & (1 << estORIRE_DTAV))
1486 update_orires_history(*fcdata.orires, &state->hist);
1489 /* ############# START The update of velocities and positions ######### */
1490 int nth = gmx_omp_nthreads_get(emntUpdate);
1492 #pragma omp parallel for num_threads(nth) schedule(static)
1493 for (int th = 0; th < nth; th++)
1497 int start_th, end_th;
1498 getThreadAtomRange(nth, th, homenr, &start_th, &end_th);
1500 const rvec* x_rvec = state->x.rvec_array();
1501 rvec* xp_rvec = xp_.rvec_array();
1502 rvec* v_rvec = state->v.rvec_array();
1503 const rvec* f_rvec = as_rvec_array(f.unpaddedConstArrayRef().data());
1505 switch (inputRecord.eI)
1508 do_update_md(start_th,
1518 inputRecord.nsttcouple,
1519 inputRecord.nstpcouple,
1523 state->nosehoover_vxi.data(),
1527 do_update_sd(start_th,
1535 inputRecord.opts.nFreeze,
1540 inputRecord.ld_seed,
1546 do_update_bd(start_th,
1554 inputRecord.opts.nFreeze,
1559 inputRecord.bd_fric,
1561 inputRecord.ld_seed,
1562 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1567 gmx_bool bExtended = (inputRecord.etc == etcNOSEHOOVER || inputRecord.epc == epcPARRINELLORAHMAN
1568 || inputRecord.epc == epcMTTK);
1570 /* assuming barostat coupled to group 0 */
1571 real alpha = 1.0 + DIM / static_cast<real>(inputRecord.opts.nrdf[0]);
1576 do_update_vv_vel(start_th,
1579 inputRecord.opts.nFreeze,
1590 do_update_vv_pos(start_th,
1593 inputRecord.opts.nFreeze,
1605 default: gmx_fatal(FARGS, "Don't know how to update coordinates");
1608 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
1612 void Update::Impl::update_for_constraint_virial(const t_inputrec& inputRecord,
1613 const t_mdatoms& md,
1614 const t_state& state,
1615 const gmx::ArrayRefWithPadding<const gmx::RVec>& f,
1616 const gmx_ekindata_t& ekind)
1618 GMX_ASSERT(inputRecord.eI == eiMD || inputRecord.eI == eiSD1,
1619 "Only leap-frog is supported here");
1621 // Cast to real for faster code, no loss in precision
1622 const real dt = inputRecord.delta_t;
1624 const int nth = gmx_omp_nthreads_get(emntUpdate);
1626 #pragma omp parallel for num_threads(nth) schedule(static)
1627 for (int th = 0; th < nth; th++)
1631 int start_th, end_th;
1632 getThreadAtomRange(nth, th, md.homenr, &start_th, &end_th);
1634 const rvec* x_rvec = state.x.rvec_array();
1635 rvec* xp_rvec = xp_.rvec_array();
1636 rvec* v_rvec = const_cast<rvec*>(state.v.rvec_array());
1637 const rvec* f_rvec = as_rvec_array(f.unpaddedConstArrayRef().data());
1639 doUpdateMDDoNotUpdateVelocities(
1640 start_th, end_th, dt, x_rvec, xp_rvec, v_rvec, f_rvec, md, ekind);
1642 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR