Revert "Eliminate mdlib/mdrun.h"
[gromacs.git] / src / gromacs / mdlib / update.cpp
blobc1f2c983a595be181b14069f3573bc68be2b4c85
1 /*
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,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 #include "gmxpre.h"
39 #include "update.h"
41 #include <cmath>
42 #include <cstdio>
44 #include <algorithm>
45 #include <memory>
47 #include "gromacs/domdec/domdec_struct.h"
48 #include "gromacs/fileio/confio.h"
49 #include "gromacs/gmxlib/network.h"
50 #include "gromacs/gmxlib/nrnb.h"
51 #include "gromacs/listed_forces/disre.h"
52 #include "gromacs/listed_forces/orires.h"
53 #include "gromacs/math/functions.h"
54 #include "gromacs/math/invertmatrix.h"
55 #include "gromacs/math/paddedvector.h"
56 #include "gromacs/math/units.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/math/vecdump.h"
59 #include "gromacs/mdlib/boxdeformation.h"
60 #include "gromacs/mdlib/constr.h"
61 #include "gromacs/mdlib/gmx_omp_nthreads.h"
62 #include "gromacs/mdlib/mdrun.h"
63 #include "gromacs/mdlib/sim_util.h"
64 #include "gromacs/mdlib/tgroup.h"
65 #include "gromacs/mdtypes/commrec.h"
66 #include "gromacs/mdtypes/group.h"
67 #include "gromacs/mdtypes/inputrec.h"
68 #include "gromacs/mdtypes/md_enums.h"
69 #include "gromacs/mdtypes/state.h"
70 #include "gromacs/pbcutil/boxutilities.h"
71 #include "gromacs/pbcutil/mshift.h"
72 #include "gromacs/pbcutil/pbc.h"
73 #include "gromacs/pulling/pull.h"
74 #include "gromacs/random/tabulatednormaldistribution.h"
75 #include "gromacs/random/threefry.h"
76 #include "gromacs/simd/simd.h"
77 #include "gromacs/timing/wallcycle.h"
78 #include "gromacs/topology/atoms.h"
79 #include "gromacs/utility/exceptions.h"
80 #include "gromacs/utility/fatalerror.h"
81 #include "gromacs/utility/futil.h"
82 #include "gromacs/utility/gmxassert.h"
83 #include "gromacs/utility/gmxomp.h"
84 #include "gromacs/utility/smalloc.h"
86 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
88 struct gmx_sd_const_t {
89 double em;
92 struct gmx_sd_sigma_t {
93 real V;
96 struct gmx_stochd_t
98 /* BD stuff */
99 std::vector<real> bd_rf;
100 /* SD stuff */
101 std::vector<gmx_sd_const_t> sdc;
102 std::vector<gmx_sd_sigma_t> sdsig;
103 /* andersen temperature control stuff */
104 std::vector<bool> randomize_group;
105 std::vector<real> boltzfac;
107 explicit gmx_stochd_t(const t_inputrec *ir);
110 //! pImpled implementation for Update
111 class Update::Impl
113 public:
114 //! Constructor
115 Impl(const t_inputrec *ir, BoxDeformation *boxDeformation);
116 //! Destructor
117 ~Impl() = default;
118 //! stochastic dynamics struct
119 std::unique_ptr<gmx_stochd_t> sd;
120 //! xprime for constraint algorithms
121 PaddedVector<RVec> xp;
122 //! Box deformation handler (or nullptr if inactive).
123 BoxDeformation *deform = nullptr;
126 Update::Update(const t_inputrec *ir, BoxDeformation *boxDeformation)
127 : impl_(new Impl(ir, boxDeformation))
130 Update::~Update()
133 gmx_stochd_t* Update::sd() const
135 return impl_->sd.get();
138 PaddedVector<RVec> * Update::xp()
140 return &impl_->xp;
143 BoxDeformation * Update::deform() const
145 return impl_->deform;
148 static bool isTemperatureCouplingStep(int64_t step, const t_inputrec *ir)
150 /* We should only couple after a step where energies were determined (for leapfrog versions)
151 or the step energies are determined, for velocity verlet versions */
152 int offset;
153 if (EI_VV(ir->eI))
155 offset = 0;
157 else
159 offset = 1;
161 return ir->etc != etcNO &&
162 (ir->nsttcouple == 1 ||
163 do_per_step(step + ir->nsttcouple - offset, ir->nsttcouple));
166 static bool isPressureCouplingStep(int64_t step, const t_inputrec *ir)
168 GMX_ASSERT(ir->epc != epcMTTK, "MTTK pressure coupling is not handled here");
170 int offset;
171 if (ir->epc == epcBERENDSEN)
173 offset = 0;
175 else
177 offset = 1;
179 /* We should only couple after a step where pressures were determined */
180 return ir->epc != etcNO &&
181 (ir->nstpcouple == 1 ||
182 do_per_step(step + ir->nstpcouple - offset, ir->nstpcouple));
185 /*! \brief Sets the velocities of virtual sites to zero */
186 static void clearVsiteVelocities(int start,
187 int nrend,
188 const unsigned short *particleType,
189 rvec * gmx_restrict v)
191 for (int a = start; a < nrend; a++)
193 if (particleType[a] == eptVSite)
195 clear_rvec(v[a]);
200 /*! \brief Sets the number of different temperature coupling values */
201 enum class NumTempScaleValues
203 single, //!< Single T-scaling value (either one group or all values =1)
204 multiple //!< Multiple T-scaling values, need to use T-group indices
207 /*! \brief Sets if to apply no or diagonal Parrinello-Rahman pressure scaling
209 * Note that this enum is only used in updateMDLeapfrogSimple(), which does
210 * not handle fully anistropic Parrinello-Rahman scaling, so we only have
211 * options \p no and \p diagonal here and no anistropic option.
213 enum class ApplyParrinelloRahmanVScaling
215 no, //!< Do not apply velocity scaling (not a PR-coupling run or step)
216 diagonal //!< Apply velocity scaling using a diagonal matrix
219 /*! \brief Integrate using leap-frog with T-scaling and optionally diagonal Parrinello-Rahman p-coupling
221 * \tparam numTempScaleValues The number of different T-couple values
222 * \tparam applyPRVScaling Apply Parrinello-Rahman velocity scaling
223 * \param[in] start Index of first atom to update
224 * \param[in] nrend Last atom to update: \p nrend - 1
225 * \param[in] dt The time step
226 * \param[in] dtPressureCouple Time step for pressure coupling
227 * \param[in] invMassPerDim 1/mass per atom and dimension
228 * \param[in] tcstat Temperature coupling information
229 * \param[in] cTC T-coupling group index per atom
230 * \param[in] pRVScaleMatrixDiagonal Parrinello-Rahman v-scale matrix diagonal
231 * \param[in] x Input coordinates
232 * \param[out] xprime Updated coordinates
233 * \param[inout] v Velocities
234 * \param[in] f Forces
236 * We expect this template to get good SIMD acceleration by most compilers,
237 * unlike the more complex general template.
238 * Note that we might get even better SIMD acceleration when we introduce
239 * aligned (and padded) memory, possibly with some hints for the compilers.
241 template<NumTempScaleValues numTempScaleValues,
242 ApplyParrinelloRahmanVScaling applyPRVScaling>
243 static void
244 updateMDLeapfrogSimple(int start,
245 int nrend,
246 real dt,
247 real dtPressureCouple,
248 const rvec * gmx_restrict invMassPerDim,
249 gmx::ArrayRef<const t_grp_tcstat> tcstat,
250 const unsigned short * cTC,
251 const rvec pRVScaleMatrixDiagonal,
252 const rvec * gmx_restrict x,
253 rvec * gmx_restrict xprime,
254 rvec * gmx_restrict v,
255 const rvec * gmx_restrict f)
257 real lambdaGroup;
259 if (numTempScaleValues == NumTempScaleValues::single)
261 lambdaGroup = tcstat[0].lambda;
264 for (int a = start; a < nrend; a++)
266 if (numTempScaleValues == NumTempScaleValues::multiple)
268 lambdaGroup = tcstat[cTC[a]].lambda;
271 for (int d = 0; d < DIM; d++)
273 /* Note that using rvec invMassPerDim results in more efficient
274 * SIMD code, but this increases the cache pressure.
275 * For large systems with PME on the CPU this slows down the
276 * (then already slow) update by 20%. If all data remains in cache,
277 * using rvec is much faster.
279 real vNew = lambdaGroup*v[a][d] + f[a][d]*invMassPerDim[a][d]*dt;
281 if (applyPRVScaling == ApplyParrinelloRahmanVScaling::diagonal)
283 vNew -= dtPressureCouple*pRVScaleMatrixDiagonal[d]*v[a][d];
285 v[a][d] = vNew;
286 xprime[a][d] = x[a][d] + vNew*dt;
291 #if GMX_SIMD && GMX_SIMD_HAVE_REAL
292 #define GMX_HAVE_SIMD_UPDATE 1
293 #else
294 #define GMX_HAVE_SIMD_UPDATE 0
295 #endif
297 #if GMX_HAVE_SIMD_UPDATE
299 /*! \brief Load (aligned) the contents of GMX_SIMD_REAL_WIDTH rvec elements sequentially into 3 SIMD registers
301 * The loaded output is:
302 * \p r0: { r[index][XX], r[index][YY], ... }
303 * \p r1: { ... }
304 * \p r2: { ..., r[index+GMX_SIMD_REAL_WIDTH-1][YY], r[index+GMX_SIMD_REAL_WIDTH-1][ZZ] }
306 * \param[in] r Real to an rvec array, has to be aligned to SIMD register width
307 * \param[in] index Index of the first rvec triplet of reals to load
308 * \param[out] r0 Pointer to first SIMD register
309 * \param[out] r1 Pointer to second SIMD register
310 * \param[out] r2 Pointer to third SIMD register
312 static inline void simdLoadRvecs(const rvec *r,
313 int index,
314 SimdReal *r0,
315 SimdReal *r1,
316 SimdReal *r2)
318 const real *realPtr = r[index];
320 GMX_ASSERT(isSimdAligned(realPtr), "Pointer should be SIMD aligned");
322 *r0 = simdLoad(realPtr + 0*GMX_SIMD_REAL_WIDTH);
323 *r1 = simdLoad(realPtr + 1*GMX_SIMD_REAL_WIDTH);
324 *r2 = simdLoad(realPtr + 2*GMX_SIMD_REAL_WIDTH);
327 /*! \brief Store (aligned) 3 SIMD registers sequentially to GMX_SIMD_REAL_WIDTH rvec elements
329 * The stored output is:
330 * \p r[index] = { { r0[0], r0[1], ... }
331 * ...
332 * \p r[index+GMX_SIMD_REAL_WIDTH-1] = { ... , r2[GMX_SIMD_REAL_WIDTH-2], r2[GMX_SIMD_REAL_WIDTH-1] }
334 * \param[out] r Pointer to an rvec array, has to be aligned to SIMD register width
335 * \param[in] index Index of the first rvec triplet of reals to store to
336 * \param[in] r0 First SIMD register
337 * \param[in] r1 Second SIMD register
338 * \param[in] r2 Third SIMD register
340 static inline void simdStoreRvecs(rvec *r,
341 int index,
342 SimdReal r0,
343 SimdReal r1,
344 SimdReal r2)
346 real *realPtr = r[index];
348 GMX_ASSERT(isSimdAligned(realPtr), "Pointer should be SIMD aligned");
350 store(realPtr + 0*GMX_SIMD_REAL_WIDTH, r0);
351 store(realPtr + 1*GMX_SIMD_REAL_WIDTH, r1);
352 store(realPtr + 2*GMX_SIMD_REAL_WIDTH, r2);
355 /*! \brief Integrate using leap-frog with single group T-scaling and SIMD
357 * \param[in] start Index of first atom to update
358 * \param[in] nrend Last atom to update: \p nrend - 1
359 * \param[in] dt The time step
360 * \param[in] invMass 1/mass per atom
361 * \param[in] tcstat Temperature coupling information
362 * \param[in] x Input coordinates
363 * \param[out] xprime Updated coordinates
364 * \param[inout] v Velocities
365 * \param[in] f Forces
367 static void
368 updateMDLeapfrogSimpleSimd(int start,
369 int nrend,
370 real dt,
371 const real * gmx_restrict invMass,
372 gmx::ArrayRef<const t_grp_tcstat> tcstat,
373 const rvec * gmx_restrict x,
374 rvec * gmx_restrict xprime,
375 rvec * gmx_restrict v,
376 const rvec * gmx_restrict f)
378 SimdReal timestep(dt);
379 SimdReal lambdaSystem(tcstat[0].lambda);
381 /* We declare variables here, since code is often slower when declaring them inside the loop */
383 /* Note: We should implement a proper PaddedVector, so we don't need this check */
384 GMX_ASSERT(isSimdAligned(invMass), "invMass should be aligned");
386 for (int a = start; a < nrend; a += GMX_SIMD_REAL_WIDTH)
388 SimdReal invMass0, invMass1, invMass2;
389 expandScalarsToTriplets(simdLoad(invMass + a),
390 &invMass0, &invMass1, &invMass2);
392 SimdReal v0, v1, v2;
393 SimdReal f0, f1, f2;
394 simdLoadRvecs(v, a, &v0, &v1, &v2);
395 simdLoadRvecs(f, a, &f0, &f1, &f2);
397 v0 = fma(f0*invMass0, timestep, lambdaSystem*v0);
398 v1 = fma(f1*invMass1, timestep, lambdaSystem*v1);
399 v2 = fma(f2*invMass2, timestep, lambdaSystem*v2);
401 simdStoreRvecs(v, a, v0, v1, v2);
403 SimdReal x0, x1, x2;
404 simdLoadRvecs(x, a, &x0, &x1, &x2);
406 SimdReal xprime0 = fma(v0, timestep, x0);
407 SimdReal xprime1 = fma(v1, timestep, x1);
408 SimdReal xprime2 = fma(v2, timestep, x2);
410 simdStoreRvecs(xprime, a, xprime0, xprime1, xprime2);
414 #endif // GMX_HAVE_SIMD_UPDATE
416 /*! \brief Sets the NEMD acceleration type */
417 enum class AccelerationType
419 none, group, cosine
422 /*! \brief Integrate using leap-frog with support for everything
424 * \tparam accelerationType Type of NEMD acceleration
425 * \param[in] start Index of first atom to update
426 * \param[in] nrend Last atom to update: \p nrend - 1
427 * \param[in] doNoseHoover If to apply Nose-Hoover T-coupling
428 * \param[in] dt The time step
429 * \param[in] dtPressureCouple Time step for pressure coupling, is 0 when no pressure coupling should be applied at this step
430 * \param[in] ir The input parameter record
431 * \param[in] md Atom properties
432 * \param[in] ekind Kinetic energy data
433 * \param[in] box The box dimensions
434 * \param[in] x Input coordinates
435 * \param[out] xprime Updated coordinates
436 * \param[inout] v Velocities
437 * \param[in] f Forces
438 * \param[in] nh_vxi Nose-Hoover velocity scaling factors
439 * \param[in] M Parrinello-Rahman scaling matrix
441 template<AccelerationType accelerationType>
442 static void
443 updateMDLeapfrogGeneral(int start,
444 int nrend,
445 bool doNoseHoover,
446 real dt,
447 real dtPressureCouple,
448 const t_inputrec * ir,
449 const t_mdatoms * md,
450 const gmx_ekindata_t * ekind,
451 const matrix box,
452 const rvec * gmx_restrict x,
453 rvec * gmx_restrict xprime,
454 rvec * gmx_restrict v,
455 const rvec * gmx_restrict f,
456 const double * gmx_restrict nh_vxi,
457 const matrix M)
459 /* This is a version of the leap-frog integrator that supports
460 * all combinations of T-coupling, P-coupling and NEMD.
461 * Nose-Hoover uses the reversible leap-frog integrator from
462 * Holian et al. Phys Rev E 52(3) : 2338, 1995
465 gmx::ArrayRef<const t_grp_tcstat> tcstat = ekind->tcstat;
466 gmx::ArrayRef<const t_grp_acc> grpstat = ekind->grpstat;
467 const unsigned short * cTC = md->cTC;
468 const unsigned short * cACC = md->cACC;
469 const rvec * accel = ir->opts.acc;
471 const rvec * gmx_restrict invMassPerDim = md->invMassPerDim;
473 /* Initialize group values, changed later when multiple groups are used */
474 int ga = 0;
475 int gt = 0;
476 real factorNH = 0;
478 real omega_Z = 2*static_cast<real>(M_PI)/box[ZZ][ZZ];
480 for (int n = start; n < nrend; n++)
482 if (cTC)
484 gt = cTC[n];
486 real lg = tcstat[gt].lambda;
488 rvec vRel;
489 real cosineZ, vCosine;
490 #ifdef __INTEL_COMPILER
491 #pragma warning( disable : 280 )
492 #endif
493 switch (accelerationType)
495 case AccelerationType::none:
496 copy_rvec(v[n], vRel);
497 break;
498 case AccelerationType::group:
499 if (cACC)
501 ga = cACC[n];
503 /* Avoid scaling the group velocity */
504 rvec_sub(v[n], grpstat[ga].u, vRel);
505 break;
506 case AccelerationType::cosine:
507 cosineZ = std::cos(x[n][ZZ]*omega_Z);
508 vCosine = cosineZ*ekind->cosacc.vcos;
509 /* Avoid scaling the cosine profile velocity */
510 copy_rvec(v[n], vRel);
511 vRel[XX] -= vCosine;
512 break;
515 if (doNoseHoover)
517 /* Here we account for multiple time stepping, by increasing
518 * the Nose-Hoover correction by nsttcouple
520 factorNH = 0.5*ir->nsttcouple*dt*nh_vxi[gt];
523 for (int d = 0; d < DIM; d++)
525 real vNew =
526 (lg*vRel[d] + (f[n][d]*invMassPerDim[n][d]*dt - factorNH*vRel[d]
527 - dtPressureCouple*iprod(M[d], vRel)))/(1 + factorNH);
528 switch (accelerationType)
530 case AccelerationType::none:
531 break;
532 case AccelerationType::group:
533 /* Add back the mean velocity and apply acceleration */
534 vNew += grpstat[ga].u[d] + accel[ga][d]*dt;
535 break;
536 case AccelerationType::cosine:
537 if (d == XX)
539 /* Add back the mean velocity and apply acceleration */
540 vNew += vCosine + cosineZ*ekind->cosacc.cos_accel*dt;
542 break;
544 v[n][d] = vNew;
545 xprime[n][d] = x[n][d] + vNew*dt;
550 /*! \brief Handles the Leap-frog MD x and v integration */
551 static void do_update_md(int start,
552 int nrend,
553 int64_t step,
554 real dt,
555 const t_inputrec * ir,
556 const t_mdatoms * md,
557 const gmx_ekindata_t * ekind,
558 const matrix box,
559 const rvec * gmx_restrict x,
560 rvec * gmx_restrict xprime,
561 rvec * gmx_restrict v,
562 const rvec * gmx_restrict f,
563 const double * gmx_restrict nh_vxi,
564 const matrix M)
566 GMX_ASSERT(nrend == start || xprime != x, "For SIMD optimization certain compilers need to have xprime != x");
568 /* Note: Berendsen pressure scaling is handled after do_update_md() */
569 bool doTempCouple = isTemperatureCouplingStep(step, ir);
570 bool doNoseHoover = (ir->etc == etcNOSEHOOVER && doTempCouple);
571 bool doParrinelloRahman = (ir->epc == epcPARRINELLORAHMAN && isPressureCouplingStep(step, ir));
572 bool doPROffDiagonal = (doParrinelloRahman && (M[YY][XX] != 0 || M[ZZ][XX] != 0 || M[ZZ][YY] != 0));
574 real dtPressureCouple = (doParrinelloRahman ? ir->nstpcouple*dt : 0);
576 /* NEMD (also cosine) acceleration is applied in updateMDLeapFrogGeneral */
577 bool doAcceleration = (ekind->bNEMD || ekind->cosacc.cos_accel != 0);
579 if (doNoseHoover || doPROffDiagonal || doAcceleration)
581 matrix stepM;
582 if (!doParrinelloRahman)
584 /* We should not apply PR scaling at this step */
585 clear_mat(stepM);
587 else
589 copy_mat(M, stepM);
592 if (!doAcceleration)
594 updateMDLeapfrogGeneral<AccelerationType::none>
595 (start, nrend, doNoseHoover, dt, dtPressureCouple,
596 ir, md, ekind, box, x, xprime, v, f, nh_vxi, stepM);
598 else if (ekind->bNEMD)
600 updateMDLeapfrogGeneral<AccelerationType::group>
601 (start, nrend, doNoseHoover, dt, dtPressureCouple,
602 ir, md, ekind, box, x, xprime, v, f, nh_vxi, stepM);
604 else
606 updateMDLeapfrogGeneral<AccelerationType::cosine>
607 (start, nrend, doNoseHoover, dt, dtPressureCouple,
608 ir, md, ekind, box, x, xprime, v, f, nh_vxi, stepM);
611 else
613 /* Use a simple and thus more efficient integration loop. */
614 /* The simple loop does not check for particle type (so it can
615 * be vectorized), which means we need to clear the velocities
616 * of virtual sites in advance, when present. Note that vsite
617 * velocities are computed after update and constraints from
618 * their displacement.
620 if (md->haveVsites)
622 /* Note: The overhead of this loop is completely neligible */
623 clearVsiteVelocities(start, nrend, md->ptype, v);
626 /* We determine if we have a single T-coupling lambda value for all
627 * atoms. That allows for better SIMD acceleration in the template.
628 * If we do not do temperature coupling (in the run or this step),
629 * all scaling values are 1, so we effectively have a single value.
631 bool haveSingleTempScaleValue = (!doTempCouple || ekind->ngtc == 1);
633 /* Extract some pointers needed by all cases */
634 const unsigned short *cTC = md->cTC;
635 gmx::ArrayRef<const t_grp_tcstat> tcstat = ekind->tcstat;
636 const rvec *invMassPerDim = md->invMassPerDim;
638 if (doParrinelloRahman)
640 GMX_ASSERT(!doPROffDiagonal, "updateMDLeapfrogSimple only support diagonal Parrinello-Rahman scaling matrices");
642 rvec diagM;
643 for (int d = 0; d < DIM; d++)
645 diagM[d] = M[d][d];
648 if (haveSingleTempScaleValue)
650 updateMDLeapfrogSimple
651 <NumTempScaleValues::single,
652 ApplyParrinelloRahmanVScaling::diagonal>
653 (start, nrend, dt, dtPressureCouple,
654 invMassPerDim, tcstat, cTC, diagM, x, xprime, v, f);
656 else
658 updateMDLeapfrogSimple
659 <NumTempScaleValues::multiple,
660 ApplyParrinelloRahmanVScaling::diagonal>
661 (start, nrend, dt, dtPressureCouple,
662 invMassPerDim, tcstat, cTC, diagM, x, xprime, v, f);
665 else
667 if (haveSingleTempScaleValue)
669 /* Note that modern compilers are pretty good at vectorizing
670 * updateMDLeapfrogSimple(). But the SIMD version will still
671 * be faster because invMass lowers the cache pressure
672 * compared to invMassPerDim.
674 #if GMX_HAVE_SIMD_UPDATE
675 /* Check if we can use invmass instead of invMassPerDim */
676 if (!md->havePartiallyFrozenAtoms)
678 updateMDLeapfrogSimpleSimd
679 (start, nrend, dt, md->invmass, tcstat, x, xprime, v, f);
681 else
682 #endif
684 updateMDLeapfrogSimple
685 <NumTempScaleValues::single,
686 ApplyParrinelloRahmanVScaling::no>
687 (start, nrend, dt, dtPressureCouple,
688 invMassPerDim, tcstat, cTC, nullptr, x, xprime, v, f);
691 else
693 updateMDLeapfrogSimple
694 <NumTempScaleValues::multiple,
695 ApplyParrinelloRahmanVScaling::no>
696 (start, nrend, dt, dtPressureCouple,
697 invMassPerDim, tcstat, cTC, nullptr, x, xprime, v, f);
703 static void do_update_vv_vel(int start, int nrend, real dt,
704 const rvec accel[], const ivec nFreeze[], const real invmass[],
705 const unsigned short ptype[], const unsigned short cFREEZE[],
706 const unsigned short cACC[], rvec v[], const rvec f[],
707 gmx_bool bExtended, real veta, real alpha)
709 int gf = 0, ga = 0;
710 int n, d;
711 real g, mv1, mv2;
713 if (bExtended)
715 g = 0.25*dt*veta*alpha;
716 mv1 = std::exp(-g);
717 mv2 = gmx::series_sinhx(g);
719 else
721 mv1 = 1.0;
722 mv2 = 1.0;
724 for (n = start; n < nrend; n++)
726 real w_dt = invmass[n]*dt;
727 if (cFREEZE)
729 gf = cFREEZE[n];
731 if (cACC)
733 ga = cACC[n];
736 for (d = 0; d < DIM; d++)
738 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
740 v[n][d] = mv1*(mv1*v[n][d] + 0.5*(w_dt*mv2*f[n][d]))+0.5*accel[ga][d]*dt;
742 else
744 v[n][d] = 0.0;
748 } /* do_update_vv_vel */
750 static void do_update_vv_pos(int start, int nrend, real dt,
751 const ivec nFreeze[],
752 const unsigned short ptype[], const unsigned short cFREEZE[],
753 const rvec x[], rvec xprime[], const rvec v[],
754 gmx_bool bExtended, real veta)
756 int gf = 0;
757 int n, d;
758 real g, mr1, mr2;
760 /* Would it make more sense if Parrinello-Rahman was put here? */
761 if (bExtended)
763 g = 0.5*dt*veta;
764 mr1 = std::exp(g);
765 mr2 = gmx::series_sinhx(g);
767 else
769 mr1 = 1.0;
770 mr2 = 1.0;
773 for (n = start; n < nrend; n++)
776 if (cFREEZE)
778 gf = cFREEZE[n];
781 for (d = 0; d < DIM; d++)
783 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
785 xprime[n][d] = mr1*(mr1*x[n][d]+mr2*dt*v[n][d]);
787 else
789 xprime[n][d] = x[n][d];
793 } /* do_update_vv_pos */
795 gmx_stochd_t::gmx_stochd_t(const t_inputrec *ir)
797 const t_grpopts *opts = &ir->opts;
798 const int ngtc = opts->ngtc;
800 if (ir->eI == eiBD)
802 bd_rf.resize(ngtc);
804 else if (EI_SD(ir->eI))
806 sdc.resize(ngtc);
807 sdsig.resize(ngtc);
809 for (int gt = 0; gt < ngtc; gt++)
811 if (opts->tau_t[gt] > 0)
813 sdc[gt].em = std::exp(-ir->delta_t/opts->tau_t[gt]);
815 else
817 /* No friction and noise on this group */
818 sdc[gt].em = 1;
822 else if (ETC_ANDERSEN(ir->etc))
824 randomize_group.resize(ngtc);
825 boltzfac.resize(ngtc);
827 /* for now, assume that all groups, if randomized, are randomized at the same rate, i.e. tau_t is the same. */
828 /* since constraint groups don't necessarily match up with temperature groups! This is checked in readir.c */
830 for (int gt = 0; gt < ngtc; gt++)
832 real reft = std::max<real>(0, opts->ref_t[gt]);
833 if ((opts->tau_t[gt] > 0) && (reft > 0)) /* tau_t or ref_t = 0 means that no randomization is done */
835 randomize_group[gt] = true;
836 boltzfac[gt] = BOLTZ*opts->ref_t[gt];
838 else
840 randomize_group[gt] = false;
846 void update_temperature_constants(gmx_stochd_t *sd, const t_inputrec *ir)
848 if (ir->eI == eiBD)
850 if (ir->bd_fric != 0)
852 for (int gt = 0; gt < ir->opts.ngtc; gt++)
854 sd->bd_rf[gt] = std::sqrt(2.0*BOLTZ*ir->opts.ref_t[gt]/(ir->bd_fric*ir->delta_t));
857 else
859 for (int gt = 0; gt < ir->opts.ngtc; gt++)
861 sd->bd_rf[gt] = std::sqrt(2.0*BOLTZ*ir->opts.ref_t[gt]);
865 if (ir->eI == eiSD1)
867 for (int gt = 0; gt < ir->opts.ngtc; gt++)
869 real kT = BOLTZ*ir->opts.ref_t[gt];
870 /* The mass is accounted for later, since this differs per atom */
871 sd->sdsig[gt].V = std::sqrt(kT*(1 - sd->sdc[gt].em * sd->sdc[gt].em));
876 Update::Impl::Impl(const t_inputrec *ir, BoxDeformation *boxDeformation)
878 sd = std::make_unique<gmx_stochd_t>(ir);
879 update_temperature_constants(sd.get(), ir);
880 xp.resizeWithPadding(0);
881 deform = boxDeformation;
884 void Update::setNumAtoms(int nAtoms)
887 impl_->xp.resizeWithPadding(nAtoms);
890 /*! \brief Sets the SD update type */
891 enum class SDUpdate : int
893 ForcesOnly, FrictionAndNoiseOnly, Combined
896 /*! \brief SD integrator update
898 * Two phases are required in the general case of a constrained
899 * update, the first phase from the contribution of forces, before
900 * applying constraints, and then a second phase applying the friction
901 * and noise, and then further constraining. For details, see
902 * Goga2012.
904 * Without constraints, the two phases can be combined, for
905 * efficiency.
907 * Thus three instantiations of this templated function will be made,
908 * two with only one contribution, and one with both contributions. */
909 template <SDUpdate updateType>
910 static void
911 doSDUpdateGeneral(const gmx_stochd_t &sd,
912 int start, int nrend, real dt,
913 const rvec accel[], const ivec nFreeze[],
914 const real invmass[], const unsigned short ptype[],
915 const unsigned short cFREEZE[], const unsigned short cACC[],
916 const unsigned short cTC[],
917 const rvec x[], rvec xprime[], rvec v[], const rvec f[],
918 int64_t step, int seed, const int *gatindex)
920 // cTC, cACC and cFreeze can be nullptr any time, but various
921 // instantiations do not make sense with particular pointer
922 // values.
923 if (updateType == SDUpdate::ForcesOnly)
925 GMX_ASSERT(f != nullptr, "SD update with only forces requires forces");
926 GMX_ASSERT(cTC == nullptr, "SD update with only forces cannot handle temperature groups");
928 if (updateType == SDUpdate::FrictionAndNoiseOnly)
930 GMX_ASSERT(f == nullptr, "SD update with only noise cannot handle forces");
931 GMX_ASSERT(cACC == nullptr, "SD update with only noise cannot handle acceleration groups");
933 if (updateType == SDUpdate::Combined)
935 GMX_ASSERT(f != nullptr, "SD update with forces and noise requires forces");
938 // Even 0 bits internal counter gives 2x64 ints (more than enough for three table lookups)
939 gmx::ThreeFry2x64<0> rng(seed, gmx::RandomDomain::UpdateCoordinates);
940 gmx::TabulatedNormalDistribution<real, 14> dist;
942 for (int n = start; n < nrend; n++)
944 int globalAtomIndex = gatindex ? gatindex[n] : n;
945 rng.restart(step, globalAtomIndex);
946 dist.reset();
948 real inverseMass = invmass[n];
949 real invsqrtMass = std::sqrt(inverseMass);
951 int freezeGroup = cFREEZE ? cFREEZE[n] : 0;
952 int accelerationGroup = cACC ? cACC[n] : 0;
953 int temperatureGroup = cTC ? cTC[n] : 0;
955 for (int d = 0; d < DIM; d++)
957 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[freezeGroup][d])
959 if (updateType == SDUpdate::ForcesOnly)
961 real vn = v[n][d] + (inverseMass*f[n][d] + accel[accelerationGroup][d])*dt;
962 v[n][d] = vn;
963 // Simple position update.
964 xprime[n][d] = x[n][d] + v[n][d]*dt;
966 else if (updateType == SDUpdate::FrictionAndNoiseOnly)
968 real vn = v[n][d];
969 v[n][d] = (vn*sd.sdc[temperatureGroup].em +
970 invsqrtMass*sd.sdsig[temperatureGroup].V*dist(rng));
971 // The previous phase already updated the
972 // positions with a full v*dt term that must
973 // now be half removed.
974 xprime[n][d] = xprime[n][d] + 0.5*(v[n][d] - vn)*dt;
976 else
978 real vn = v[n][d] + (inverseMass*f[n][d] + accel[accelerationGroup][d])*dt;
979 v[n][d] = (vn*sd.sdc[temperatureGroup].em +
980 invsqrtMass*sd.sdsig[temperatureGroup].V*dist(rng));
981 // Here we include half of the friction+noise
982 // update of v into the position update.
983 xprime[n][d] = x[n][d] + 0.5*(vn + v[n][d])*dt;
986 else
988 // When using constraints, the update is split into
989 // two phases, but we only need to zero the update of
990 // virtual, shell or frozen particles in at most one
991 // of the phases.
992 if (updateType != SDUpdate::FrictionAndNoiseOnly)
994 v[n][d] = 0.0;
995 xprime[n][d] = x[n][d];
1002 static void do_update_bd(int start, int nrend, real dt,
1003 const ivec nFreeze[],
1004 const real invmass[], const unsigned short ptype[],
1005 const unsigned short cFREEZE[], const unsigned short cTC[],
1006 const rvec x[], rvec xprime[], rvec v[],
1007 const rvec f[], real friction_coefficient,
1008 const real *rf, int64_t step, int seed,
1009 const int* gatindex)
1011 /* note -- these appear to be full step velocities . . . */
1012 int gf = 0, gt = 0;
1013 real vn;
1014 real invfr = 0;
1015 int n, d;
1016 // Use 1 bit of internal counters to give us 2*2 64-bits values per stream
1017 // Each 64-bit value is enough for 4 normal distribution table numbers.
1018 gmx::ThreeFry2x64<0> rng(seed, gmx::RandomDomain::UpdateCoordinates);
1019 gmx::TabulatedNormalDistribution<real, 14> dist;
1021 if (friction_coefficient != 0)
1023 invfr = 1.0/friction_coefficient;
1026 for (n = start; (n < nrend); n++)
1028 int ng = gatindex ? gatindex[n] : n;
1030 rng.restart(step, ng);
1031 dist.reset();
1033 if (cFREEZE)
1035 gf = cFREEZE[n];
1037 if (cTC)
1039 gt = cTC[n];
1041 for (d = 0; (d < DIM); d++)
1043 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
1045 if (friction_coefficient != 0)
1047 vn = invfr*f[n][d] + rf[gt]*dist(rng);
1049 else
1051 /* NOTE: invmass = 2/(mass*friction_constant*dt) */
1052 vn = 0.5*invmass[n]*f[n][d]*dt
1053 + std::sqrt(0.5*invmass[n])*rf[gt]*dist(rng);
1056 v[n][d] = vn;
1057 xprime[n][d] = x[n][d]+vn*dt;
1059 else
1061 v[n][d] = 0.0;
1062 xprime[n][d] = x[n][d];
1068 static void calc_ke_part_normal(const rvec v[], const t_grpopts *opts, const t_mdatoms *md,
1069 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel)
1071 int g;
1072 gmx::ArrayRef<t_grp_tcstat> tcstat = ekind->tcstat;
1073 gmx::ArrayRef<t_grp_acc> grpstat = ekind->grpstat;
1074 int nthread, thread;
1076 /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin. Leap with AveVel is also
1077 an option, but not supported now.
1078 bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
1081 /* group velocities are calculated in update_ekindata and
1082 * accumulated in acumulate_groups.
1083 * Now the partial global and groups ekin.
1085 for (g = 0; (g < opts->ngtc); g++)
1087 copy_mat(tcstat[g].ekinh, tcstat[g].ekinh_old);
1088 if (bEkinAveVel)
1090 clear_mat(tcstat[g].ekinf);
1091 tcstat[g].ekinscalef_nhc = 1.0; /* need to clear this -- logic is complicated! */
1093 else
1095 clear_mat(tcstat[g].ekinh);
1098 ekind->dekindl_old = ekind->dekindl;
1099 nthread = gmx_omp_nthreads_get(emntUpdate);
1101 #pragma omp parallel for num_threads(nthread) schedule(static)
1102 for (thread = 0; thread < nthread; thread++)
1104 // This OpenMP only loops over arrays and does not call any functions
1105 // or memory allocation. It should not be able to throw, so for now
1106 // we do not need a try/catch wrapper.
1107 int start_t, end_t, n;
1108 int ga, gt;
1109 rvec v_corrt;
1110 real hm;
1111 int d, m;
1112 matrix *ekin_sum;
1113 real *dekindl_sum;
1115 start_t = ((thread+0)*md->homenr)/nthread;
1116 end_t = ((thread+1)*md->homenr)/nthread;
1118 ekin_sum = ekind->ekin_work[thread];
1119 dekindl_sum = ekind->dekindl_work[thread];
1121 for (gt = 0; gt < opts->ngtc; gt++)
1123 clear_mat(ekin_sum[gt]);
1125 *dekindl_sum = 0.0;
1127 ga = 0;
1128 gt = 0;
1129 for (n = start_t; n < end_t; n++)
1131 if (md->cACC)
1133 ga = md->cACC[n];
1135 if (md->cTC)
1137 gt = md->cTC[n];
1139 hm = 0.5*md->massT[n];
1141 for (d = 0; (d < DIM); d++)
1143 v_corrt[d] = v[n][d] - grpstat[ga].u[d];
1145 for (d = 0; (d < DIM); d++)
1147 for (m = 0; (m < DIM); m++)
1149 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1150 ekin_sum[gt][m][d] += hm*v_corrt[m]*v_corrt[d];
1153 if (md->nMassPerturbed && md->bPerturbed[n])
1155 *dekindl_sum +=
1156 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1161 ekind->dekindl = 0;
1162 for (thread = 0; thread < nthread; thread++)
1164 for (g = 0; g < opts->ngtc; g++)
1166 if (bEkinAveVel)
1168 m_add(tcstat[g].ekinf, ekind->ekin_work[thread][g],
1169 tcstat[g].ekinf);
1171 else
1173 m_add(tcstat[g].ekinh, ekind->ekin_work[thread][g],
1174 tcstat[g].ekinh);
1178 ekind->dekindl += *ekind->dekindl_work[thread];
1181 inc_nrnb(nrnb, eNR_EKIN, md->homenr);
1184 static void calc_ke_part_visc(const matrix box, const rvec x[], const rvec v[],
1185 const t_grpopts *opts, const t_mdatoms *md,
1186 gmx_ekindata_t *ekind,
1187 t_nrnb *nrnb, gmx_bool bEkinAveVel)
1189 int start = 0, homenr = md->homenr;
1190 int g, d, n, m, gt = 0;
1191 rvec v_corrt;
1192 real hm;
1193 gmx::ArrayRef<t_grp_tcstat> tcstat = ekind->tcstat;
1194 t_cos_acc *cosacc = &(ekind->cosacc);
1195 real dekindl;
1196 real fac, cosz;
1197 double mvcos;
1199 for (g = 0; g < opts->ngtc; g++)
1201 copy_mat(ekind->tcstat[g].ekinh, ekind->tcstat[g].ekinh_old);
1202 clear_mat(ekind->tcstat[g].ekinh);
1204 ekind->dekindl_old = ekind->dekindl;
1206 fac = 2*M_PI/box[ZZ][ZZ];
1207 mvcos = 0;
1208 dekindl = 0;
1209 for (n = start; n < start+homenr; n++)
1211 if (md->cTC)
1213 gt = md->cTC[n];
1215 hm = 0.5*md->massT[n];
1217 /* Note that the times of x and v differ by half a step */
1218 /* MRS -- would have to be changed for VV */
1219 cosz = std::cos(fac*x[n][ZZ]);
1220 /* Calculate the amplitude of the new velocity profile */
1221 mvcos += 2*cosz*md->massT[n]*v[n][XX];
1223 copy_rvec(v[n], v_corrt);
1224 /* Subtract the profile for the kinetic energy */
1225 v_corrt[XX] -= cosz*cosacc->vcos;
1226 for (d = 0; (d < DIM); d++)
1228 for (m = 0; (m < DIM); m++)
1230 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1231 if (bEkinAveVel)
1233 tcstat[gt].ekinf[m][d] += hm*v_corrt[m]*v_corrt[d];
1235 else
1237 tcstat[gt].ekinh[m][d] += hm*v_corrt[m]*v_corrt[d];
1241 if (md->nPerturbed && md->bPerturbed[n])
1243 /* The minus sign here might be confusing.
1244 * The kinetic contribution from dH/dl doesn't come from
1245 * d m(l)/2 v^2 / dl, but rather from d p^2/2m(l) / dl,
1246 * where p are the momenta. The difference is only a minus sign.
1248 dekindl -= 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1251 ekind->dekindl = dekindl;
1252 cosacc->mvcos = mvcos;
1254 inc_nrnb(nrnb, eNR_EKIN, homenr);
1257 void calc_ke_part(const t_state *state, const t_grpopts *opts, const t_mdatoms *md,
1258 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel)
1260 if (ekind->cosacc.cos_accel == 0)
1262 calc_ke_part_normal(state->v.rvec_array(), opts, md, ekind, nrnb, bEkinAveVel);
1264 else
1266 calc_ke_part_visc(state->box, state->x.rvec_array(), state->v.rvec_array(), opts, md, ekind, nrnb, bEkinAveVel);
1270 extern void init_ekinstate(ekinstate_t *ekinstate, const t_inputrec *ir)
1272 ekinstate->ekin_n = ir->opts.ngtc;
1273 snew(ekinstate->ekinh, ekinstate->ekin_n);
1274 snew(ekinstate->ekinf, ekinstate->ekin_n);
1275 snew(ekinstate->ekinh_old, ekinstate->ekin_n);
1276 ekinstate->ekinscalef_nhc.resize(ekinstate->ekin_n);
1277 ekinstate->ekinscaleh_nhc.resize(ekinstate->ekin_n);
1278 ekinstate->vscale_nhc.resize(ekinstate->ekin_n);
1279 ekinstate->dekindl = 0;
1280 ekinstate->mvcos = 0;
1283 void update_ekinstate(ekinstate_t *ekinstate, const gmx_ekindata_t *ekind)
1285 int i;
1287 for (i = 0; i < ekinstate->ekin_n; i++)
1289 copy_mat(ekind->tcstat[i].ekinh, ekinstate->ekinh[i]);
1290 copy_mat(ekind->tcstat[i].ekinf, ekinstate->ekinf[i]);
1291 copy_mat(ekind->tcstat[i].ekinh_old, ekinstate->ekinh_old[i]);
1292 ekinstate->ekinscalef_nhc[i] = ekind->tcstat[i].ekinscalef_nhc;
1293 ekinstate->ekinscaleh_nhc[i] = ekind->tcstat[i].ekinscaleh_nhc;
1294 ekinstate->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
1297 copy_mat(ekind->ekin, ekinstate->ekin_total);
1298 ekinstate->dekindl = ekind->dekindl;
1299 ekinstate->mvcos = ekind->cosacc.mvcos;
1303 void restore_ekinstate_from_state(const t_commrec *cr,
1304 gmx_ekindata_t *ekind, const ekinstate_t *ekinstate)
1306 int i, n;
1308 if (MASTER(cr))
1310 for (i = 0; i < ekinstate->ekin_n; i++)
1312 copy_mat(ekinstate->ekinh[i], ekind->tcstat[i].ekinh);
1313 copy_mat(ekinstate->ekinf[i], ekind->tcstat[i].ekinf);
1314 copy_mat(ekinstate->ekinh_old[i], ekind->tcstat[i].ekinh_old);
1315 ekind->tcstat[i].ekinscalef_nhc = ekinstate->ekinscalef_nhc[i];
1316 ekind->tcstat[i].ekinscaleh_nhc = ekinstate->ekinscaleh_nhc[i];
1317 ekind->tcstat[i].vscale_nhc = ekinstate->vscale_nhc[i];
1320 copy_mat(ekinstate->ekin_total, ekind->ekin);
1322 ekind->dekindl = ekinstate->dekindl;
1323 ekind->cosacc.mvcos = ekinstate->mvcos;
1324 n = ekinstate->ekin_n;
1327 if (PAR(cr))
1329 gmx_bcast(sizeof(n), &n, cr);
1330 for (i = 0; i < n; i++)
1332 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh[0][0]),
1333 ekind->tcstat[i].ekinh[0], cr);
1334 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinf[0][0]),
1335 ekind->tcstat[i].ekinf[0], cr);
1336 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh_old[0][0]),
1337 ekind->tcstat[i].ekinh_old[0], cr);
1339 gmx_bcast(sizeof(ekind->tcstat[i].ekinscalef_nhc),
1340 &(ekind->tcstat[i].ekinscalef_nhc), cr);
1341 gmx_bcast(sizeof(ekind->tcstat[i].ekinscaleh_nhc),
1342 &(ekind->tcstat[i].ekinscaleh_nhc), cr);
1343 gmx_bcast(sizeof(ekind->tcstat[i].vscale_nhc),
1344 &(ekind->tcstat[i].vscale_nhc), cr);
1346 gmx_bcast(DIM*DIM*sizeof(ekind->ekin[0][0]),
1347 ekind->ekin[0], cr);
1349 gmx_bcast(sizeof(ekind->dekindl), &ekind->dekindl, cr);
1350 gmx_bcast(sizeof(ekind->cosacc.mvcos), &ekind->cosacc.mvcos, cr);
1354 void update_tcouple(int64_t step,
1355 const t_inputrec *inputrec,
1356 t_state *state,
1357 gmx_ekindata_t *ekind,
1358 const t_extmass *MassQ,
1359 const t_mdatoms *md)
1362 bool doTemperatureCoupling = false;
1364 /* if using vv with trotter decomposition methods, we do this elsewhere in the code */
1365 if (!(EI_VV(inputrec->eI) &&
1366 (inputrecNvtTrotter(inputrec) || inputrecNptTrotter(inputrec) || inputrecNphTrotter(inputrec))))
1368 doTemperatureCoupling = isTemperatureCouplingStep(step, inputrec);
1371 if (doTemperatureCoupling)
1373 real dttc = inputrec->nsttcouple*inputrec->delta_t;
1375 switch (inputrec->etc)
1377 case etcNO:
1378 break;
1379 case etcBERENDSEN:
1380 berendsen_tcoupl(inputrec, ekind, dttc, state->therm_integral);
1381 break;
1382 case etcNOSEHOOVER:
1383 nosehoover_tcoupl(&(inputrec->opts), ekind, dttc,
1384 state->nosehoover_xi.data(), state->nosehoover_vxi.data(), MassQ);
1385 break;
1386 case etcVRESCALE:
1387 vrescale_tcoupl(inputrec, step, ekind, dttc,
1388 state->therm_integral.data());
1389 break;
1391 /* rescale in place here */
1392 if (EI_VV(inputrec->eI))
1394 rescale_velocities(ekind, md, 0, md->homenr, state->v.rvec_array());
1397 else
1399 /* Set the T scaling lambda to 1 to have no scaling */
1400 for (int i = 0; (i < inputrec->opts.ngtc); i++)
1402 ekind->tcstat[i].lambda = 1.0;
1407 /*! \brief Computes the atom range for a thread to operate on, ensured SIMD aligned ranges
1409 * \param[in] numThreads The number of threads to divide atoms over
1410 * \param[in] threadIndex The thread to get the range for
1411 * \param[in] numAtoms The total number of atoms (on this rank)
1412 * \param[out] startAtom The start of the atom range
1413 * \param[out] endAtom The end of the atom range, note that this is in general not a multiple of the SIMD width
1415 static void getThreadAtomRange(int numThreads, int threadIndex, int numAtoms,
1416 int *startAtom, int *endAtom)
1418 #if GMX_HAVE_SIMD_UPDATE
1419 constexpr int blockSize = GMX_SIMD_REAL_WIDTH;
1420 #else
1421 constexpr int blockSize = 1;
1422 #endif
1423 int numBlocks = (numAtoms + blockSize - 1)/blockSize;
1425 *startAtom = ((numBlocks* threadIndex )/numThreads)*blockSize;
1426 *endAtom = ((numBlocks*(threadIndex + 1))/numThreads)*blockSize;
1427 if (threadIndex == numThreads - 1)
1429 *endAtom = numAtoms;
1433 void update_pcouple_before_coordinates(FILE *fplog,
1434 int64_t step,
1435 const t_inputrec *inputrec,
1436 t_state *state,
1437 matrix parrinellorahmanMu,
1438 matrix M,
1439 gmx_bool bInitStep)
1441 /* Berendsen P-coupling is completely handled after the coordinate update.
1442 * Trotter P-coupling is handled by separate calls to trotter_update().
1444 if (inputrec->epc == epcPARRINELLORAHMAN &&
1445 isPressureCouplingStep(step, inputrec))
1447 real dtpc = inputrec->nstpcouple*inputrec->delta_t;
1449 parrinellorahman_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev,
1450 state->box, state->box_rel, state->boxv,
1451 M, parrinellorahmanMu, bInitStep);
1455 void constrain_velocities(int64_t step,
1456 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1457 t_state *state,
1458 tensor vir_part,
1459 gmx::Constraints *constr,
1460 gmx_bool bCalcVir,
1461 bool do_log,
1462 bool do_ene)
1464 if (!constr)
1466 return;
1470 * Steps (7C, 8C)
1471 * APPLY CONSTRAINTS:
1472 * BLOCK SHAKE
1476 tensor vir_con;
1478 /* clear out constraints before applying */
1479 clear_mat(vir_part);
1481 /* Constrain the coordinates upd->xp */
1482 constr->apply(do_log, do_ene,
1483 step, 1, 1.0,
1484 state->x.rvec_array(), state->v.rvec_array(), state->v.rvec_array(),
1485 state->box,
1486 state->lambda[efptBONDED], dvdlambda,
1487 nullptr, bCalcVir ? &vir_con : nullptr, ConstraintVariable::Velocities);
1489 if (bCalcVir)
1491 m_add(vir_part, vir_con, vir_part);
1496 void constrain_coordinates(int64_t step,
1497 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1498 t_state *state,
1499 tensor vir_part,
1500 Update *upd,
1501 gmx::Constraints *constr,
1502 gmx_bool bCalcVir,
1503 bool do_log,
1504 bool do_ene)
1506 if (!constr)
1508 return;
1512 tensor vir_con;
1514 /* clear out constraints before applying */
1515 clear_mat(vir_part);
1517 /* Constrain the coordinates upd->xp */
1518 constr->apply(do_log, do_ene,
1519 step, 1, 1.0,
1520 state->x.rvec_array(), upd->xp()->rvec_array(), nullptr,
1521 state->box,
1522 state->lambda[efptBONDED], dvdlambda,
1523 as_rvec_array(state->v.data()), bCalcVir ? &vir_con : nullptr, ConstraintVariable::Positions);
1525 if (bCalcVir)
1527 m_add(vir_part, vir_con, vir_part);
1532 void
1533 update_sd_second_half(int64_t step,
1534 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1535 const t_inputrec *inputrec, /* input record and box stuff */
1536 const t_mdatoms *md,
1537 t_state *state,
1538 const t_commrec *cr,
1539 t_nrnb *nrnb,
1540 gmx_wallcycle_t wcycle,
1541 Update *upd,
1542 gmx::Constraints *constr,
1543 bool do_log,
1544 bool do_ene)
1546 if (!constr)
1548 return;
1550 if (inputrec->eI == eiSD1)
1552 int nth, th;
1553 int homenr = md->homenr;
1555 /* Cast delta_t from double to real to make the integrators faster.
1556 * The only reason for having delta_t double is to get accurate values
1557 * for t=delta_t*step when step is larger than float precision.
1558 * For integration dt the accuracy of real suffices, since with
1559 * integral += dt*integrand the increment is nearly always (much) smaller
1560 * than the integral (and the integrand has real precision).
1562 real dt = inputrec->delta_t;
1564 wallcycle_start(wcycle, ewcUPDATE);
1566 nth = gmx_omp_nthreads_get(emntUpdate);
1568 #pragma omp parallel for num_threads(nth) schedule(static)
1569 for (th = 0; th < nth; th++)
1573 int start_th, end_th;
1574 getThreadAtomRange(nth, th, homenr, &start_th, &end_th);
1576 doSDUpdateGeneral<SDUpdate::FrictionAndNoiseOnly>
1577 (*upd->sd(),
1578 start_th, end_th, dt,
1579 inputrec->opts.acc, inputrec->opts.nFreeze,
1580 md->invmass, md->ptype,
1581 md->cFREEZE, nullptr, md->cTC,
1582 state->x.rvec_array(), upd->xp()->rvec_array(),
1583 state->v.rvec_array(), nullptr,
1584 step, inputrec->ld_seed,
1585 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1587 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1589 inc_nrnb(nrnb, eNR_UPDATE, homenr);
1590 wallcycle_stop(wcycle, ewcUPDATE);
1592 /* Constrain the coordinates upd->xp for half a time step */
1593 constr->apply(do_log, do_ene,
1594 step, 1, 0.5,
1595 state->x.rvec_array(), upd->xp()->rvec_array(), nullptr,
1596 state->box,
1597 state->lambda[efptBONDED], dvdlambda,
1598 as_rvec_array(state->v.data()), nullptr, ConstraintVariable::Positions);
1602 void finish_update(const t_inputrec *inputrec, /* input record and box stuff */
1603 const t_mdatoms *md,
1604 t_state *state,
1605 const t_graph *graph,
1606 t_nrnb *nrnb,
1607 gmx_wallcycle_t wcycle,
1608 Update *upd,
1609 const gmx::Constraints *constr)
1611 int homenr = md->homenr;
1613 /* We must always unshift after updating coordinates; if we did not shake
1614 x was shifted in do_force */
1616 /* NOTE Currently we always integrate to a temporary buffer and
1617 * then copy the results back. */
1619 wallcycle_start_nocount(wcycle, ewcUPDATE);
1621 if (md->cFREEZE != nullptr && constr != nullptr)
1623 /* If we have atoms that are frozen along some, but not all
1624 * dimensions, then any constraints will have moved them also along
1625 * the frozen dimensions. To freeze such degrees of freedom
1626 * we copy them back here to later copy them forward. It would
1627 * be more elegant and slightly more efficient to copies zero
1628 * times instead of twice, but the graph case below prevents this.
1630 const ivec *nFreeze = inputrec->opts.nFreeze;
1631 bool partialFreezeAndConstraints = false;
1632 for (int g = 0; g < inputrec->opts.ngfrz; g++)
1634 int numFreezeDim = nFreeze[g][XX] + nFreeze[g][YY] + nFreeze[g][ZZ];
1635 if (numFreezeDim > 0 && numFreezeDim < 3)
1637 partialFreezeAndConstraints = true;
1640 if (partialFreezeAndConstraints)
1642 auto xp = makeArrayRef(*upd->xp()).subArray(0, homenr);
1643 auto x = makeConstArrayRef(state->x).subArray(0, homenr);
1644 for (int i = 0; i < homenr; i++)
1646 int g = md->cFREEZE[i];
1648 for (int d = 0; d < DIM; d++)
1650 if (nFreeze[g][d])
1652 xp[i][d] = x[i][d];
1659 if (graph && (graph->nnodes > 0))
1661 unshift_x(graph, state->box, state->x.rvec_array(), upd->xp()->rvec_array());
1662 if (TRICLINIC(state->box))
1664 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
1666 else
1668 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
1671 else
1673 auto xp = makeConstArrayRef(*upd->xp()).subArray(0, homenr);
1674 auto x = makeArrayRef(state->x).subArray(0, homenr);
1675 #ifndef __clang_analyzer__
1676 int gmx_unused nth = gmx_omp_nthreads_get(emntUpdate);
1677 #endif
1678 #pragma omp parallel for num_threads(nth) schedule(static)
1679 for (int i = 0; i < homenr; i++)
1681 // Trivial statement, does not throw
1682 x[i] = xp[i];
1685 wallcycle_stop(wcycle, ewcUPDATE);
1687 /* ############# END the update of velocities and positions ######### */
1690 void update_pcouple_after_coordinates(FILE *fplog,
1691 int64_t step,
1692 const t_inputrec *inputrec,
1693 const t_mdatoms *md,
1694 const matrix pressure,
1695 const matrix forceVirial,
1696 const matrix constraintVirial,
1697 const matrix parrinellorahmanMu,
1698 t_state *state,
1699 t_nrnb *nrnb,
1700 Update *upd)
1702 int start = 0;
1703 int homenr = md->homenr;
1705 /* Cast to real for faster code, no loss in precision (see comment above) */
1706 real dt = inputrec->delta_t;
1709 /* now update boxes */
1710 switch (inputrec->epc)
1712 case (epcNO):
1713 break;
1714 case (epcBERENDSEN):
1715 if (isPressureCouplingStep(step, inputrec))
1717 real dtpc = inputrec->nstpcouple*dt;
1718 matrix mu;
1719 berendsen_pcoupl(fplog, step, inputrec, dtpc,
1720 pressure, state->box,
1721 forceVirial, constraintVirial,
1722 mu, &state->baros_integral);
1723 berendsen_pscale(inputrec, mu, state->box, state->box_rel,
1724 start, homenr, state->x.rvec_array(),
1725 md->cFREEZE, nrnb);
1727 break;
1728 case (epcPARRINELLORAHMAN):
1729 if (isPressureCouplingStep(step, inputrec))
1731 /* The box velocities were updated in do_pr_pcoupl,
1732 * but we dont change the box vectors until we get here
1733 * since we need to be able to shift/unshift above.
1735 real dtpc = inputrec->nstpcouple*dt;
1736 for (int i = 0; i < DIM; i++)
1738 for (int m = 0; m <= i; m++)
1740 state->box[i][m] += dtpc*state->boxv[i][m];
1743 preserve_box_shape(inputrec, state->box_rel, state->box);
1745 /* Scale the coordinates */
1746 auto x = state->x.rvec_array();
1747 for (int n = start; n < start + homenr; n++)
1749 tmvmul_ur0(parrinellorahmanMu, x[n], x[n]);
1752 break;
1753 case (epcMTTK):
1754 switch (inputrec->epct)
1756 case (epctISOTROPIC):
1757 /* DIM * eta = ln V. so DIM*eta_new = DIM*eta_old + DIM*dt*veta =>
1758 ln V_new = ln V_old + 3*dt*veta => V_new = V_old*exp(3*dt*veta) =>
1759 Side length scales as exp(veta*dt) */
1761 msmul(state->box, std::exp(state->veta*dt), state->box);
1763 /* Relate veta to boxv. veta = d(eta)/dT = (1/DIM)*1/V dV/dT.
1764 o If we assume isotropic scaling, and box length scaling
1765 factor L, then V = L^DIM (det(M)). So dV/dt = DIM
1766 L^(DIM-1) dL/dt det(M), and veta = (1/L) dL/dt. The
1767 determinant of B is L^DIM det(M), and the determinant
1768 of dB/dt is (dL/dT)^DIM det (M). veta will be
1769 (det(dB/dT)/det(B))^(1/3). Then since M =
1770 B_new*(vol_new)^(1/3), dB/dT_new = (veta_new)*B(new). */
1772 msmul(state->box, state->veta, state->boxv);
1773 break;
1774 default:
1775 break;
1777 break;
1778 default:
1779 break;
1782 if (upd->deform())
1784 auto localX = makeArrayRef(state->x).subArray(start, homenr);
1785 upd->deform()->apply(localX, state->box, step);
1789 void update_coords(int64_t step,
1790 const t_inputrec *inputrec, /* input record and box stuff */
1791 const t_mdatoms *md,
1792 t_state *state,
1793 gmx::ArrayRefWithPadding<gmx::RVec> f,
1794 const t_fcdata *fcd,
1795 const gmx_ekindata_t *ekind,
1796 const matrix M,
1797 Update *upd,
1798 int UpdatePart,
1799 const t_commrec *cr, /* these shouldn't be here -- need to think about it */
1800 const gmx::Constraints *constr)
1802 gmx_bool bDoConstr = (nullptr != constr);
1804 /* Running the velocity half does nothing except for velocity verlet */
1805 if ((UpdatePart == etrtVELOCITY1 || UpdatePart == etrtVELOCITY2) &&
1806 !EI_VV(inputrec->eI))
1808 gmx_incons("update_coords called for velocity without VV integrator");
1811 int homenr = md->homenr;
1813 /* Cast to real for faster code, no loss in precision (see comment above) */
1814 real dt = inputrec->delta_t;
1816 /* We need to update the NMR restraint history when time averaging is used */
1817 if (state->flags & (1<<estDISRE_RM3TAV))
1819 update_disres_history(fcd, &state->hist);
1821 if (state->flags & (1<<estORIRE_DTAV))
1823 update_orires_history(fcd, &state->hist);
1826 /* ############# START The update of velocities and positions ######### */
1827 int nth = gmx_omp_nthreads_get(emntUpdate);
1829 #pragma omp parallel for num_threads(nth) schedule(static)
1830 for (int th = 0; th < nth; th++)
1834 int start_th, end_th;
1835 getThreadAtomRange(nth, th, homenr, &start_th, &end_th);
1837 const rvec *x_rvec = state->x.rvec_array();
1838 rvec *xp_rvec = upd->xp()->rvec_array();
1839 rvec *v_rvec = state->v.rvec_array();
1840 const rvec *f_rvec = as_rvec_array(f.unpaddedArrayRef().data());
1842 switch (inputrec->eI)
1844 case (eiMD):
1845 do_update_md(start_th, end_th, step, dt,
1846 inputrec, md, ekind, state->box,
1847 x_rvec, xp_rvec, v_rvec, f_rvec,
1848 state->nosehoover_vxi.data(), M);
1849 break;
1850 case (eiSD1):
1851 if (bDoConstr)
1853 // With constraints, the SD update is done in 2 parts
1854 doSDUpdateGeneral<SDUpdate::ForcesOnly>
1855 (*upd->sd(),
1856 start_th, end_th, dt,
1857 inputrec->opts.acc, inputrec->opts.nFreeze,
1858 md->invmass, md->ptype,
1859 md->cFREEZE, md->cACC, nullptr,
1860 x_rvec, xp_rvec, v_rvec, f_rvec,
1861 step, inputrec->ld_seed, nullptr);
1863 else
1865 doSDUpdateGeneral<SDUpdate::Combined>
1866 (*upd->sd(),
1867 start_th, end_th, dt,
1868 inputrec->opts.acc, inputrec->opts.nFreeze,
1869 md->invmass, md->ptype,
1870 md->cFREEZE, md->cACC, md->cTC,
1871 x_rvec, xp_rvec, v_rvec, f_rvec,
1872 step, inputrec->ld_seed,
1873 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1875 break;
1876 case (eiBD):
1877 do_update_bd(start_th, end_th, dt,
1878 inputrec->opts.nFreeze, md->invmass, md->ptype,
1879 md->cFREEZE, md->cTC,
1880 x_rvec, xp_rvec, v_rvec, f_rvec,
1881 inputrec->bd_fric,
1882 upd->sd()->bd_rf.data(),
1883 step, inputrec->ld_seed, DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1884 break;
1885 case (eiVV):
1886 case (eiVVAK):
1888 gmx_bool bExtended = (inputrec->etc == etcNOSEHOOVER ||
1889 inputrec->epc == epcPARRINELLORAHMAN ||
1890 inputrec->epc == epcMTTK);
1892 /* assuming barostat coupled to group 0 */
1893 real alpha = 1.0 + DIM/static_cast<real>(inputrec->opts.nrdf[0]);
1894 switch (UpdatePart)
1896 case etrtVELOCITY1:
1897 case etrtVELOCITY2:
1898 do_update_vv_vel(start_th, end_th, dt,
1899 inputrec->opts.acc, inputrec->opts.nFreeze,
1900 md->invmass, md->ptype,
1901 md->cFREEZE, md->cACC,
1902 v_rvec, f_rvec,
1903 bExtended, state->veta, alpha);
1904 break;
1905 case etrtPOSITION:
1906 do_update_vv_pos(start_th, end_th, dt,
1907 inputrec->opts.nFreeze,
1908 md->ptype, md->cFREEZE,
1909 x_rvec, xp_rvec, v_rvec,
1910 bExtended, state->veta);
1911 break;
1913 break;
1915 default:
1916 gmx_fatal(FARGS, "Don't know how to update coordinates");
1919 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1924 extern gmx_bool update_randomize_velocities(const t_inputrec *ir, int64_t step, const t_commrec *cr,
1925 const t_mdatoms *md,
1926 gmx::ArrayRef<gmx::RVec> v,
1927 const Update *upd,
1928 const gmx::Constraints *constr)
1931 real rate = (ir->delta_t)/ir->opts.tau_t[0];
1933 if (ir->etc == etcANDERSEN && constr != nullptr)
1935 /* Currently, Andersen thermostat does not support constrained
1936 systems. Functionality exists in the andersen_tcoupl
1937 function in GROMACS 4.5.7 to allow this combination. That
1938 code could be ported to the current random-number
1939 generation approach, but has not yet been done because of
1940 lack of time and resources. */
1941 gmx_fatal(FARGS, "Normal Andersen is currently not supported with constraints, use massive Andersen instead");
1944 /* proceed with andersen if 1) it's fixed probability per
1945 particle andersen or 2) it's massive andersen and it's tau_t/dt */
1946 if ((ir->etc == etcANDERSEN) || do_per_step(step, roundToInt(1.0/rate)))
1948 andersen_tcoupl(ir, step, cr, md, v, rate,
1949 upd->sd()->randomize_group,
1950 upd->sd()->boltzfac);
1951 return TRUE;
1953 return FALSE;