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.
46 #include "gromacs/compat/make_unique.h"
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
99 std::vector<real> bd_rf;
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 gmx_stochd_t(const t_inputrec *ir);
112 std::unique_ptr<gmx_stochd_t> sd;
113 /* xprime for constraint algorithms */
114 PaddedVector<gmx::RVec> xp;
116 /* Variables for the deform algorithm */
117 int64_t deformref_step;
118 matrix deformref_box;
120 //! Box deformation handler (or nullptr if inactive).
121 gmx::BoxDeformation *deform;
124 static bool isTemperatureCouplingStep(int64_t step, const t_inputrec *ir)
126 /* We should only couple after a step where energies were determined (for leapfrog versions)
127 or the step energies are determined, for velocity verlet versions */
137 return ir->etc != etcNO &&
138 (ir->nsttcouple == 1 ||
139 do_per_step(step + ir->nsttcouple - offset, ir->nsttcouple));
142 static bool isPressureCouplingStep(int64_t step, const t_inputrec *ir)
144 GMX_ASSERT(ir->epc != epcMTTK, "MTTK pressure coupling is not handled here");
147 if (ir->epc == epcBERENDSEN)
155 /* We should only couple after a step where pressures were determined */
156 return ir->epc != etcNO &&
157 (ir->nstpcouple == 1 ||
158 do_per_step(step + ir->nstpcouple - offset, ir->nstpcouple));
161 /*! \brief Sets the velocities of virtual sites to zero */
162 static void clearVsiteVelocities(int start,
164 const unsigned short *particleType,
165 rvec * gmx_restrict v)
167 for (int a = start; a < nrend; a++)
169 if (particleType[a] == eptVSite)
176 /*! \brief Sets the number of different temperature coupling values */
177 enum class NumTempScaleValues
179 single, //!< Single T-scaling value (either one group or all values =1)
180 multiple //!< Multiple T-scaling values, need to use T-group indices
183 /*! \brief Sets if to apply no or diagonal Parrinello-Rahman pressure scaling
185 * Note that this enum is only used in updateMDLeapfrogSimple(), which does
186 * not handle fully anistropic Parrinello-Rahman scaling, so we only have
187 * options \p no and \p diagonal here and no anistropic option.
189 enum class ApplyParrinelloRahmanVScaling
191 no, //!< Do not apply velocity scaling (not a PR-coupling run or step)
192 diagonal //!< Apply velocity scaling using a diagonal matrix
195 /*! \brief Integrate using leap-frog with T-scaling and optionally diagonal Parrinello-Rahman p-coupling
197 * \tparam numTempScaleValues The number of different T-couple values
198 * \tparam applyPRVScaling Apply Parrinello-Rahman velocity scaling
199 * \param[in] start Index of first atom to update
200 * \param[in] nrend Last atom to update: \p nrend - 1
201 * \param[in] dt The time step
202 * \param[in] dtPressureCouple Time step for pressure coupling
203 * \param[in] invMassPerDim 1/mass per atom and dimension
204 * \param[in] tcstat Temperature coupling information
205 * \param[in] cTC T-coupling group index per atom
206 * \param[in] pRVScaleMatrixDiagonal Parrinello-Rahman v-scale matrix diagonal
207 * \param[in] x Input coordinates
208 * \param[out] xprime Updated coordinates
209 * \param[inout] v Velocities
210 * \param[in] f Forces
212 * We expect this template to get good SIMD acceleration by most compilers,
213 * unlike the more complex general template.
214 * Note that we might get even better SIMD acceleration when we introduce
215 * aligned (and padded) memory, possibly with some hints for the compilers.
217 template<NumTempScaleValues numTempScaleValues,
218 ApplyParrinelloRahmanVScaling applyPRVScaling>
220 updateMDLeapfrogSimple(int start,
223 real dtPressureCouple,
224 const rvec * gmx_restrict invMassPerDim,
225 const t_grp_tcstat * tcstat,
226 const unsigned short * cTC,
227 const rvec pRVScaleMatrixDiagonal,
228 const rvec * gmx_restrict x,
229 rvec * gmx_restrict xprime,
230 rvec * gmx_restrict v,
231 const rvec * gmx_restrict f)
235 if (numTempScaleValues == NumTempScaleValues::single)
237 lambdaGroup = tcstat[0].lambda;
240 for (int a = start; a < nrend; a++)
242 if (numTempScaleValues == NumTempScaleValues::multiple)
244 lambdaGroup = tcstat[cTC[a]].lambda;
247 for (int d = 0; d < DIM; d++)
249 /* Note that using rvec invMassPerDim results in more efficient
250 * SIMD code, but this increases the cache pressure.
251 * For large systems with PME on the CPU this slows down the
252 * (then already slow) update by 20%. If all data remains in cache,
253 * using rvec is much faster.
255 real vNew = lambdaGroup*v[a][d] + f[a][d]*invMassPerDim[a][d]*dt;
257 if (applyPRVScaling == ApplyParrinelloRahmanVScaling::diagonal)
259 vNew -= dtPressureCouple*pRVScaleMatrixDiagonal[d]*v[a][d];
262 xprime[a][d] = x[a][d] + vNew*dt;
267 #if GMX_SIMD && GMX_SIMD_HAVE_REAL
268 #define GMX_HAVE_SIMD_UPDATE 1
270 #define GMX_HAVE_SIMD_UPDATE 0
273 #if GMX_HAVE_SIMD_UPDATE
275 /*! \brief Load (aligned) the contents of GMX_SIMD_REAL_WIDTH rvec elements sequentially into 3 SIMD registers
277 * The loaded output is:
278 * \p r0: { r[index][XX], r[index][YY], ... }
280 * \p r2: { ..., r[index+GMX_SIMD_REAL_WIDTH-1][YY], r[index+GMX_SIMD_REAL_WIDTH-1][ZZ] }
282 * \param[in] r Real to an rvec array, has to be aligned to SIMD register width
283 * \param[in] index Index of the first rvec triplet of reals to load
284 * \param[out] r0 Pointer to first SIMD register
285 * \param[out] r1 Pointer to second SIMD register
286 * \param[out] r2 Pointer to third SIMD register
288 static inline void simdLoadRvecs(const rvec *r,
294 const real *realPtr = r[index];
296 GMX_ASSERT(isSimdAligned(realPtr), "Pointer should be SIMD aligned");
298 *r0 = simdLoad(realPtr + 0*GMX_SIMD_REAL_WIDTH);
299 *r1 = simdLoad(realPtr + 1*GMX_SIMD_REAL_WIDTH);
300 *r2 = simdLoad(realPtr + 2*GMX_SIMD_REAL_WIDTH);
303 /*! \brief Store (aligned) 3 SIMD registers sequentially to GMX_SIMD_REAL_WIDTH rvec elements
305 * The stored output is:
306 * \p r[index] = { { r0[0], r0[1], ... }
308 * \p r[index+GMX_SIMD_REAL_WIDTH-1] = { ... , r2[GMX_SIMD_REAL_WIDTH-2], r2[GMX_SIMD_REAL_WIDTH-1] }
310 * \param[out] r Pointer to an rvec array, has to be aligned to SIMD register width
311 * \param[in] index Index of the first rvec triplet of reals to store to
312 * \param[in] r0 First SIMD register
313 * \param[in] r1 Second SIMD register
314 * \param[in] r2 Third SIMD register
316 static inline void simdStoreRvecs(rvec *r,
322 real *realPtr = r[index];
324 GMX_ASSERT(isSimdAligned(realPtr), "Pointer should be SIMD aligned");
326 store(realPtr + 0*GMX_SIMD_REAL_WIDTH, r0);
327 store(realPtr + 1*GMX_SIMD_REAL_WIDTH, r1);
328 store(realPtr + 2*GMX_SIMD_REAL_WIDTH, r2);
331 /*! \brief Integrate using leap-frog with single group T-scaling and SIMD
333 * \param[in] start Index of first atom to update
334 * \param[in] nrend Last atom to update: \p nrend - 1
335 * \param[in] dt The time step
336 * \param[in] invMass 1/mass per atom
337 * \param[in] tcstat Temperature coupling information
338 * \param[in] x Input coordinates
339 * \param[out] xprime Updated coordinates
340 * \param[inout] v Velocities
341 * \param[in] f Forces
344 updateMDLeapfrogSimpleSimd(int start,
347 const real * gmx_restrict invMass,
348 const t_grp_tcstat * tcstat,
349 const rvec * gmx_restrict x,
350 rvec * gmx_restrict xprime,
351 rvec * gmx_restrict v,
352 const rvec * gmx_restrict f)
354 SimdReal timestep(dt);
355 SimdReal lambdaSystem(tcstat[0].lambda);
357 /* We declare variables here, since code is often slower when declaring them inside the loop */
359 /* Note: We should implement a proper PaddedVector, so we don't need this check */
360 GMX_ASSERT(isSimdAligned(invMass), "invMass should be aligned");
362 for (int a = start; a < nrend; a += GMX_SIMD_REAL_WIDTH)
364 SimdReal invMass0, invMass1, invMass2;
365 expandScalarsToTriplets(simdLoad(invMass + a),
366 &invMass0, &invMass1, &invMass2);
370 simdLoadRvecs(v, a, &v0, &v1, &v2);
371 simdLoadRvecs(f, a, &f0, &f1, &f2);
373 v0 = fma(f0*invMass0, timestep, lambdaSystem*v0);
374 v1 = fma(f1*invMass1, timestep, lambdaSystem*v1);
375 v2 = fma(f2*invMass2, timestep, lambdaSystem*v2);
377 simdStoreRvecs(v, a, v0, v1, v2);
380 simdLoadRvecs(x, a, &x0, &x1, &x2);
382 SimdReal xprime0 = fma(v0, timestep, x0);
383 SimdReal xprime1 = fma(v1, timestep, x1);
384 SimdReal xprime2 = fma(v2, timestep, x2);
386 simdStoreRvecs(xprime, a, xprime0, xprime1, xprime2);
390 #endif // GMX_HAVE_SIMD_UPDATE
392 /*! \brief Sets the NEMD acceleration type */
393 enum class AccelerationType
398 /*! \brief Integrate using leap-frog with support for everything
400 * \tparam accelerationType Type of NEMD acceleration
401 * \param[in] start Index of first atom to update
402 * \param[in] nrend Last atom to update: \p nrend - 1
403 * \param[in] doNoseHoover If to apply Nose-Hoover T-coupling
404 * \param[in] dt The time step
405 * \param[in] dtPressureCouple Time step for pressure coupling, is 0 when no pressure coupling should be applied at this step
406 * \param[in] ir The input parameter record
407 * \param[in] md Atom properties
408 * \param[in] ekind Kinetic energy data
409 * \param[in] box The box dimensions
410 * \param[in] x Input coordinates
411 * \param[out] xprime Updated coordinates
412 * \param[inout] v Velocities
413 * \param[in] f Forces
414 * \param[in] nh_vxi Nose-Hoover velocity scaling factors
415 * \param[in] M Parrinello-Rahman scaling matrix
417 template<AccelerationType accelerationType>
419 updateMDLeapfrogGeneral(int start,
423 real dtPressureCouple,
424 const t_inputrec * ir,
425 const t_mdatoms * md,
426 const gmx_ekindata_t * ekind,
428 const rvec * gmx_restrict x,
429 rvec * gmx_restrict xprime,
430 rvec * gmx_restrict v,
431 const rvec * gmx_restrict f,
432 const double * gmx_restrict nh_vxi,
435 /* This is a version of the leap-frog integrator that supports
436 * all combinations of T-coupling, P-coupling and NEMD.
437 * Nose-Hoover uses the reversible leap-frog integrator from
438 * Holian et al. Phys Rev E 52(3) : 2338, 1995
441 const unsigned short * cTC = md->cTC;
442 const t_grp_tcstat * tcstat = ekind->tcstat;
444 const unsigned short * cACC = md->cACC;
445 const rvec * accel = ir->opts.acc;
446 const t_grp_acc * grpstat = ekind->grpstat;
448 const rvec * gmx_restrict invMassPerDim = md->invMassPerDim;
450 /* Initialize group values, changed later when multiple groups are used */
455 real omega_Z = 2*static_cast<real>(M_PI)/box[ZZ][ZZ];
457 for (int n = start; n < nrend; n++)
463 real lg = tcstat[gt].lambda;
466 real cosineZ, vCosine;
467 #ifdef __INTEL_COMPILER
468 #pragma warning( disable : 280 )
470 switch (accelerationType)
472 case AccelerationType::none:
473 copy_rvec(v[n], vRel);
475 case AccelerationType::group:
480 /* Avoid scaling the group velocity */
481 rvec_sub(v[n], grpstat[ga].u, vRel);
483 case AccelerationType::cosine:
484 cosineZ = std::cos(x[n][ZZ]*omega_Z);
485 vCosine = cosineZ*ekind->cosacc.vcos;
486 /* Avoid scaling the cosine profile velocity */
487 copy_rvec(v[n], vRel);
494 /* Here we account for multiple time stepping, by increasing
495 * the Nose-Hoover correction by nsttcouple
497 factorNH = 0.5*ir->nsttcouple*dt*nh_vxi[gt];
500 for (int d = 0; d < DIM; d++)
503 (lg*vRel[d] + (f[n][d]*invMassPerDim[n][d]*dt - factorNH*vRel[d]
504 - dtPressureCouple*iprod(M[d], vRel)))/(1 + factorNH);
505 switch (accelerationType)
507 case AccelerationType::none:
509 case AccelerationType::group:
510 /* Add back the mean velocity and apply acceleration */
511 vNew += grpstat[ga].u[d] + accel[ga][d]*dt;
513 case AccelerationType::cosine:
516 /* Add back the mean velocity and apply acceleration */
517 vNew += vCosine + cosineZ*ekind->cosacc.cos_accel*dt;
522 xprime[n][d] = x[n][d] + vNew*dt;
527 /*! \brief Handles the Leap-frog MD x and v integration */
528 static void do_update_md(int start,
532 const t_inputrec * ir,
533 const t_mdatoms * md,
534 const gmx_ekindata_t * ekind,
536 const rvec * gmx_restrict x,
537 rvec * gmx_restrict xprime,
538 rvec * gmx_restrict v,
539 const rvec * gmx_restrict f,
540 const double * gmx_restrict nh_vxi,
543 GMX_ASSERT(nrend == start || xprime != x, "For SIMD optimization certain compilers need to have xprime != x");
545 /* Note: Berendsen pressure scaling is handled after do_update_md() */
546 bool doTempCouple = isTemperatureCouplingStep(step, ir);
547 bool doNoseHoover = (ir->etc == etcNOSEHOOVER && doTempCouple);
548 bool doParrinelloRahman = (ir->epc == epcPARRINELLORAHMAN && isPressureCouplingStep(step, ir));
549 bool doPROffDiagonal = (doParrinelloRahman && (M[YY][XX] != 0 || M[ZZ][XX] != 0 || M[ZZ][YY] != 0));
551 real dtPressureCouple = (doParrinelloRahman ? ir->nstpcouple*dt : 0);
553 /* NEMD (also cosine) acceleration is applied in updateMDLeapFrogGeneral */
554 bool doAcceleration = (ekind->bNEMD || ekind->cosacc.cos_accel != 0);
556 if (doNoseHoover || doPROffDiagonal || doAcceleration)
559 if (!doParrinelloRahman)
561 /* We should not apply PR scaling at this step */
571 updateMDLeapfrogGeneral<AccelerationType::none>
572 (start, nrend, doNoseHoover, dt, dtPressureCouple,
573 ir, md, ekind, box, x, xprime, v, f, nh_vxi, stepM);
575 else if (ekind->bNEMD)
577 updateMDLeapfrogGeneral<AccelerationType::group>
578 (start, nrend, doNoseHoover, dt, dtPressureCouple,
579 ir, md, ekind, box, x, xprime, v, f, nh_vxi, stepM);
583 updateMDLeapfrogGeneral<AccelerationType::cosine>
584 (start, nrend, doNoseHoover, dt, dtPressureCouple,
585 ir, md, ekind, box, x, xprime, v, f, nh_vxi, stepM);
590 /* Use a simple and thus more efficient integration loop. */
591 /* The simple loop does not check for particle type (so it can
592 * be vectorized), which means we need to clear the velocities
593 * of virtual sites in advance, when present. Note that vsite
594 * velocities are computed after update and constraints from
595 * their displacement.
599 /* Note: The overhead of this loop is completely neligible */
600 clearVsiteVelocities(start, nrend, md->ptype, v);
603 /* We determine if we have a single T-coupling lambda value for all
604 * atoms. That allows for better SIMD acceleration in the template.
605 * If we do not do temperature coupling (in the run or this step),
606 * all scaling values are 1, so we effectively have a single value.
608 bool haveSingleTempScaleValue = (!doTempCouple || ekind->ngtc == 1);
610 /* Extract some pointers needed by all cases */
611 const unsigned short *cTC = md->cTC;
612 const t_grp_tcstat *tcstat = ekind->tcstat;
613 const rvec *invMassPerDim = md->invMassPerDim;
615 if (doParrinelloRahman)
617 GMX_ASSERT(!doPROffDiagonal, "updateMDLeapfrogSimple only support diagonal Parrinello-Rahman scaling matrices");
620 for (int d = 0; d < DIM; d++)
625 if (haveSingleTempScaleValue)
627 updateMDLeapfrogSimple
628 <NumTempScaleValues::single,
629 ApplyParrinelloRahmanVScaling::diagonal>
630 (start, nrend, dt, dtPressureCouple,
631 invMassPerDim, tcstat, cTC, diagM, x, xprime, v, f);
635 updateMDLeapfrogSimple
636 <NumTempScaleValues::multiple,
637 ApplyParrinelloRahmanVScaling::diagonal>
638 (start, nrend, dt, dtPressureCouple,
639 invMassPerDim, tcstat, cTC, diagM, x, xprime, v, f);
644 if (haveSingleTempScaleValue)
646 /* Note that modern compilers are pretty good at vectorizing
647 * updateMDLeapfrogSimple(). But the SIMD version will still
648 * be faster because invMass lowers the cache pressure
649 * compared to invMassPerDim.
651 #if GMX_HAVE_SIMD_UPDATE
652 /* Check if we can use invmass instead of invMassPerDim */
653 if (!md->havePartiallyFrozenAtoms)
655 updateMDLeapfrogSimpleSimd
656 (start, nrend, dt, md->invmass, tcstat, x, xprime, v, f);
661 updateMDLeapfrogSimple
662 <NumTempScaleValues::single,
663 ApplyParrinelloRahmanVScaling::no>
664 (start, nrend, dt, dtPressureCouple,
665 invMassPerDim, tcstat, cTC, nullptr, x, xprime, v, f);
670 updateMDLeapfrogSimple
671 <NumTempScaleValues::multiple,
672 ApplyParrinelloRahmanVScaling::no>
673 (start, nrend, dt, dtPressureCouple,
674 invMassPerDim, tcstat, cTC, nullptr, x, xprime, v, f);
680 static void do_update_vv_vel(int start, int nrend, real dt,
681 const rvec accel[], const ivec nFreeze[], const real invmass[],
682 const unsigned short ptype[], const unsigned short cFREEZE[],
683 const unsigned short cACC[], rvec v[], const rvec f[],
684 gmx_bool bExtended, real veta, real alpha)
692 g = 0.25*dt*veta*alpha;
694 mv2 = gmx::series_sinhx(g);
701 for (n = start; n < nrend; n++)
703 real w_dt = invmass[n]*dt;
713 for (d = 0; d < DIM; d++)
715 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
717 v[n][d] = mv1*(mv1*v[n][d] + 0.5*(w_dt*mv2*f[n][d]))+0.5*accel[ga][d]*dt;
725 } /* do_update_vv_vel */
727 static void do_update_vv_pos(int start, int nrend, real dt,
728 const ivec nFreeze[],
729 const unsigned short ptype[], const unsigned short cFREEZE[],
730 const rvec x[], rvec xprime[], const rvec v[],
731 gmx_bool bExtended, real veta)
737 /* Would it make more sense if Parrinello-Rahman was put here? */
742 mr2 = gmx::series_sinhx(g);
750 for (n = start; n < nrend; n++)
758 for (d = 0; d < DIM; d++)
760 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
762 xprime[n][d] = mr1*(mr1*x[n][d]+mr2*dt*v[n][d]);
766 xprime[n][d] = x[n][d];
770 } /* do_update_vv_pos */
772 gmx_stochd_t::gmx_stochd_t(const t_inputrec *ir)
774 const t_grpopts *opts = &ir->opts;
775 const int ngtc = opts->ngtc;
781 else if (EI_SD(ir->eI))
786 for (int gt = 0; gt < ngtc; gt++)
788 if (opts->tau_t[gt] > 0)
790 sdc[gt].em = std::exp(-ir->delta_t/opts->tau_t[gt]);
794 /* No friction and noise on this group */
799 else if (ETC_ANDERSEN(ir->etc))
801 randomize_group.resize(ngtc);
802 boltzfac.resize(ngtc);
804 /* for now, assume that all groups, if randomized, are randomized at the same rate, i.e. tau_t is the same. */
805 /* since constraint groups don't necessarily match up with temperature groups! This is checked in readir.c */
807 for (int gt = 0; gt < ngtc; gt++)
809 real reft = std::max<real>(0, opts->ref_t[gt]);
810 if ((opts->tau_t[gt] > 0) && (reft > 0)) /* tau_t or ref_t = 0 means that no randomization is done */
812 randomize_group[gt] = true;
813 boltzfac[gt] = BOLTZ*opts->ref_t[gt];
817 randomize_group[gt] = false;
823 void update_temperature_constants(gmx_update_t *upd, const t_inputrec *ir)
827 if (ir->bd_fric != 0)
829 for (int gt = 0; gt < ir->opts.ngtc; gt++)
831 upd->sd->bd_rf[gt] = std::sqrt(2.0*BOLTZ*ir->opts.ref_t[gt]/(ir->bd_fric*ir->delta_t));
836 for (int gt = 0; gt < ir->opts.ngtc; gt++)
838 upd->sd->bd_rf[gt] = std::sqrt(2.0*BOLTZ*ir->opts.ref_t[gt]);
844 for (int gt = 0; gt < ir->opts.ngtc; gt++)
846 real kT = BOLTZ*ir->opts.ref_t[gt];
847 /* The mass is accounted for later, since this differs per atom */
848 upd->sd->sdsig[gt].V = std::sqrt(kT*(1 - upd->sd->sdc[gt].em*upd->sd->sdc[gt].em));
853 gmx_update_t *init_update(const t_inputrec *ir,
854 gmx::BoxDeformation *deform)
856 gmx_update_t *upd = new(gmx_update_t);
858 upd->sd = gmx::compat::make_unique<gmx_stochd_t>(ir);
860 update_temperature_constants(upd, ir);
862 upd->xp.resizeWithPadding(0);
864 upd->deform = deform;
869 void update_realloc(gmx_update_t *upd, int natoms)
871 GMX_ASSERT(upd, "upd must be allocated before its fields can be reallocated");
873 upd->xp.resizeWithPadding(natoms);
876 /*! \brief Sets the SD update type */
877 enum class SDUpdate : int
879 ForcesOnly, FrictionAndNoiseOnly, Combined
882 /*! \brief SD integrator update
884 * Two phases are required in the general case of a constrained
885 * update, the first phase from the contribution of forces, before
886 * applying constraints, and then a second phase applying the friction
887 * and noise, and then further constraining. For details, see
890 * Without constraints, the two phases can be combined, for
893 * Thus three instantiations of this templated function will be made,
894 * two with only one contribution, and one with both contributions. */
895 template <SDUpdate updateType>
897 doSDUpdateGeneral(const gmx_stochd_t &sd,
898 int start, int nrend, real dt,
899 const rvec accel[], const ivec nFreeze[],
900 const real invmass[], const unsigned short ptype[],
901 const unsigned short cFREEZE[], const unsigned short cACC[],
902 const unsigned short cTC[],
903 const rvec x[], rvec xprime[], rvec v[], const rvec f[],
904 int64_t step, int seed, const int *gatindex)
906 // cTC, cACC and cFreeze can be nullptr any time, but various
907 // instantiations do not make sense with particular pointer
909 if (updateType == SDUpdate::ForcesOnly)
911 GMX_ASSERT(f != nullptr, "SD update with only forces requires forces");
912 GMX_ASSERT(cTC == nullptr, "SD update with only forces cannot handle temperature groups");
914 if (updateType == SDUpdate::FrictionAndNoiseOnly)
916 GMX_ASSERT(f == nullptr, "SD update with only noise cannot handle forces");
917 GMX_ASSERT(cACC == nullptr, "SD update with only noise cannot handle acceleration groups");
919 if (updateType == SDUpdate::Combined)
921 GMX_ASSERT(f != nullptr, "SD update with forces and noise requires forces");
924 // Even 0 bits internal counter gives 2x64 ints (more than enough for three table lookups)
925 gmx::ThreeFry2x64<0> rng(seed, gmx::RandomDomain::UpdateCoordinates);
926 gmx::TabulatedNormalDistribution<real, 14> dist;
928 for (int n = start; n < nrend; n++)
930 int globalAtomIndex = gatindex ? gatindex[n] : n;
931 rng.restart(step, globalAtomIndex);
934 real inverseMass = invmass[n];
935 real invsqrtMass = std::sqrt(inverseMass);
937 int freezeGroup = cFREEZE ? cFREEZE[n] : 0;
938 int accelerationGroup = cACC ? cACC[n] : 0;
939 int temperatureGroup = cTC ? cTC[n] : 0;
941 for (int d = 0; d < DIM; d++)
943 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[freezeGroup][d])
945 if (updateType == SDUpdate::ForcesOnly)
947 real vn = v[n][d] + (inverseMass*f[n][d] + accel[accelerationGroup][d])*dt;
949 // Simple position update.
950 xprime[n][d] = x[n][d] + v[n][d]*dt;
952 else if (updateType == SDUpdate::FrictionAndNoiseOnly)
955 v[n][d] = (vn*sd.sdc[temperatureGroup].em +
956 invsqrtMass*sd.sdsig[temperatureGroup].V*dist(rng));
957 // The previous phase already updated the
958 // positions with a full v*dt term that must
959 // now be half removed.
960 xprime[n][d] = xprime[n][d] + 0.5*(v[n][d] - vn)*dt;
964 real vn = v[n][d] + (inverseMass*f[n][d] + accel[accelerationGroup][d])*dt;
965 v[n][d] = (vn*sd.sdc[temperatureGroup].em +
966 invsqrtMass*sd.sdsig[temperatureGroup].V*dist(rng));
967 // Here we include half of the friction+noise
968 // update of v into the position update.
969 xprime[n][d] = x[n][d] + 0.5*(vn + v[n][d])*dt;
974 // When using constraints, the update is split into
975 // two phases, but we only need to zero the update of
976 // virtual, shell or frozen particles in at most one
978 if (updateType != SDUpdate::FrictionAndNoiseOnly)
981 xprime[n][d] = x[n][d];
988 static void do_update_bd(int start, int nrend, real dt,
989 const ivec nFreeze[],
990 const real invmass[], const unsigned short ptype[],
991 const unsigned short cFREEZE[], const unsigned short cTC[],
992 const rvec x[], rvec xprime[], rvec v[],
993 const rvec f[], real friction_coefficient,
994 const real *rf, int64_t step, int seed,
997 /* note -- these appear to be full step velocities . . . */
1002 // Use 1 bit of internal counters to give us 2*2 64-bits values per stream
1003 // Each 64-bit value is enough for 4 normal distribution table numbers.
1004 gmx::ThreeFry2x64<0> rng(seed, gmx::RandomDomain::UpdateCoordinates);
1005 gmx::TabulatedNormalDistribution<real, 14> dist;
1007 if (friction_coefficient != 0)
1009 invfr = 1.0/friction_coefficient;
1012 for (n = start; (n < nrend); n++)
1014 int ng = gatindex ? gatindex[n] : n;
1016 rng.restart(step, ng);
1027 for (d = 0; (d < DIM); d++)
1029 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
1031 if (friction_coefficient != 0)
1033 vn = invfr*f[n][d] + rf[gt]*dist(rng);
1037 /* NOTE: invmass = 2/(mass*friction_constant*dt) */
1038 vn = 0.5*invmass[n]*f[n][d]*dt
1039 + std::sqrt(0.5*invmass[n])*rf[gt]*dist(rng);
1043 xprime[n][d] = x[n][d]+vn*dt;
1048 xprime[n][d] = x[n][d];
1054 static void calc_ke_part_normal(const rvec v[], const t_grpopts *opts, const t_mdatoms *md,
1055 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel)
1058 t_grp_tcstat *tcstat = ekind->tcstat;
1059 t_grp_acc *grpstat = ekind->grpstat;
1060 int nthread, thread;
1062 /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin. Leap with AveVel is also
1063 an option, but not supported now.
1064 bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
1067 /* group velocities are calculated in update_ekindata and
1068 * accumulated in acumulate_groups.
1069 * Now the partial global and groups ekin.
1071 for (g = 0; (g < opts->ngtc); g++)
1073 copy_mat(tcstat[g].ekinh, tcstat[g].ekinh_old);
1076 clear_mat(tcstat[g].ekinf);
1077 tcstat[g].ekinscalef_nhc = 1.0; /* need to clear this -- logic is complicated! */
1081 clear_mat(tcstat[g].ekinh);
1084 ekind->dekindl_old = ekind->dekindl;
1085 nthread = gmx_omp_nthreads_get(emntUpdate);
1087 #pragma omp parallel for num_threads(nthread) schedule(static)
1088 for (thread = 0; thread < nthread; thread++)
1090 // This OpenMP only loops over arrays and does not call any functions
1091 // or memory allocation. It should not be able to throw, so for now
1092 // we do not need a try/catch wrapper.
1093 int start_t, end_t, n;
1101 start_t = ((thread+0)*md->homenr)/nthread;
1102 end_t = ((thread+1)*md->homenr)/nthread;
1104 ekin_sum = ekind->ekin_work[thread];
1105 dekindl_sum = ekind->dekindl_work[thread];
1107 for (gt = 0; gt < opts->ngtc; gt++)
1109 clear_mat(ekin_sum[gt]);
1115 for (n = start_t; n < end_t; n++)
1125 hm = 0.5*md->massT[n];
1127 for (d = 0; (d < DIM); d++)
1129 v_corrt[d] = v[n][d] - grpstat[ga].u[d];
1131 for (d = 0; (d < DIM); d++)
1133 for (m = 0; (m < DIM); m++)
1135 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1136 ekin_sum[gt][m][d] += hm*v_corrt[m]*v_corrt[d];
1139 if (md->nMassPerturbed && md->bPerturbed[n])
1142 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1148 for (thread = 0; thread < nthread; thread++)
1150 for (g = 0; g < opts->ngtc; g++)
1154 m_add(tcstat[g].ekinf, ekind->ekin_work[thread][g],
1159 m_add(tcstat[g].ekinh, ekind->ekin_work[thread][g],
1164 ekind->dekindl += *ekind->dekindl_work[thread];
1167 inc_nrnb(nrnb, eNR_EKIN, md->homenr);
1170 static void calc_ke_part_visc(const matrix box, const rvec x[], const rvec v[],
1171 const t_grpopts *opts, const t_mdatoms *md,
1172 gmx_ekindata_t *ekind,
1173 t_nrnb *nrnb, gmx_bool bEkinAveVel)
1175 int start = 0, homenr = md->homenr;
1176 int g, d, n, m, gt = 0;
1179 t_grp_tcstat *tcstat = ekind->tcstat;
1180 t_cos_acc *cosacc = &(ekind->cosacc);
1185 for (g = 0; g < opts->ngtc; g++)
1187 copy_mat(ekind->tcstat[g].ekinh, ekind->tcstat[g].ekinh_old);
1188 clear_mat(ekind->tcstat[g].ekinh);
1190 ekind->dekindl_old = ekind->dekindl;
1192 fac = 2*M_PI/box[ZZ][ZZ];
1195 for (n = start; n < start+homenr; n++)
1201 hm = 0.5*md->massT[n];
1203 /* Note that the times of x and v differ by half a step */
1204 /* MRS -- would have to be changed for VV */
1205 cosz = std::cos(fac*x[n][ZZ]);
1206 /* Calculate the amplitude of the new velocity profile */
1207 mvcos += 2*cosz*md->massT[n]*v[n][XX];
1209 copy_rvec(v[n], v_corrt);
1210 /* Subtract the profile for the kinetic energy */
1211 v_corrt[XX] -= cosz*cosacc->vcos;
1212 for (d = 0; (d < DIM); d++)
1214 for (m = 0; (m < DIM); m++)
1216 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1219 tcstat[gt].ekinf[m][d] += hm*v_corrt[m]*v_corrt[d];
1223 tcstat[gt].ekinh[m][d] += hm*v_corrt[m]*v_corrt[d];
1227 if (md->nPerturbed && md->bPerturbed[n])
1229 /* The minus sign here might be confusing.
1230 * The kinetic contribution from dH/dl doesn't come from
1231 * d m(l)/2 v^2 / dl, but rather from d p^2/2m(l) / dl,
1232 * where p are the momenta. The difference is only a minus sign.
1234 dekindl -= 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1237 ekind->dekindl = dekindl;
1238 cosacc->mvcos = mvcos;
1240 inc_nrnb(nrnb, eNR_EKIN, homenr);
1243 void calc_ke_part(const t_state *state, const t_grpopts *opts, const t_mdatoms *md,
1244 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel)
1246 if (ekind->cosacc.cos_accel == 0)
1248 calc_ke_part_normal(state->v.rvec_array(), opts, md, ekind, nrnb, bEkinAveVel);
1252 calc_ke_part_visc(state->box, state->x.rvec_array(), state->v.rvec_array(), opts, md, ekind, nrnb, bEkinAveVel);
1256 extern void init_ekinstate(ekinstate_t *ekinstate, const t_inputrec *ir)
1258 ekinstate->ekin_n = ir->opts.ngtc;
1259 snew(ekinstate->ekinh, ekinstate->ekin_n);
1260 snew(ekinstate->ekinf, ekinstate->ekin_n);
1261 snew(ekinstate->ekinh_old, ekinstate->ekin_n);
1262 ekinstate->ekinscalef_nhc.resize(ekinstate->ekin_n);
1263 ekinstate->ekinscaleh_nhc.resize(ekinstate->ekin_n);
1264 ekinstate->vscale_nhc.resize(ekinstate->ekin_n);
1265 ekinstate->dekindl = 0;
1266 ekinstate->mvcos = 0;
1269 void update_ekinstate(ekinstate_t *ekinstate, const gmx_ekindata_t *ekind)
1273 for (i = 0; i < ekinstate->ekin_n; i++)
1275 copy_mat(ekind->tcstat[i].ekinh, ekinstate->ekinh[i]);
1276 copy_mat(ekind->tcstat[i].ekinf, ekinstate->ekinf[i]);
1277 copy_mat(ekind->tcstat[i].ekinh_old, ekinstate->ekinh_old[i]);
1278 ekinstate->ekinscalef_nhc[i] = ekind->tcstat[i].ekinscalef_nhc;
1279 ekinstate->ekinscaleh_nhc[i] = ekind->tcstat[i].ekinscaleh_nhc;
1280 ekinstate->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
1283 copy_mat(ekind->ekin, ekinstate->ekin_total);
1284 ekinstate->dekindl = ekind->dekindl;
1285 ekinstate->mvcos = ekind->cosacc.mvcos;
1289 void restore_ekinstate_from_state(const t_commrec *cr,
1290 gmx_ekindata_t *ekind, const ekinstate_t *ekinstate)
1296 for (i = 0; i < ekinstate->ekin_n; i++)
1298 copy_mat(ekinstate->ekinh[i], ekind->tcstat[i].ekinh);
1299 copy_mat(ekinstate->ekinf[i], ekind->tcstat[i].ekinf);
1300 copy_mat(ekinstate->ekinh_old[i], ekind->tcstat[i].ekinh_old);
1301 ekind->tcstat[i].ekinscalef_nhc = ekinstate->ekinscalef_nhc[i];
1302 ekind->tcstat[i].ekinscaleh_nhc = ekinstate->ekinscaleh_nhc[i];
1303 ekind->tcstat[i].vscale_nhc = ekinstate->vscale_nhc[i];
1306 copy_mat(ekinstate->ekin_total, ekind->ekin);
1308 ekind->dekindl = ekinstate->dekindl;
1309 ekind->cosacc.mvcos = ekinstate->mvcos;
1310 n = ekinstate->ekin_n;
1315 gmx_bcast(sizeof(n), &n, cr);
1316 for (i = 0; i < n; i++)
1318 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh[0][0]),
1319 ekind->tcstat[i].ekinh[0], cr);
1320 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinf[0][0]),
1321 ekind->tcstat[i].ekinf[0], cr);
1322 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh_old[0][0]),
1323 ekind->tcstat[i].ekinh_old[0], cr);
1325 gmx_bcast(sizeof(ekind->tcstat[i].ekinscalef_nhc),
1326 &(ekind->tcstat[i].ekinscalef_nhc), cr);
1327 gmx_bcast(sizeof(ekind->tcstat[i].ekinscaleh_nhc),
1328 &(ekind->tcstat[i].ekinscaleh_nhc), cr);
1329 gmx_bcast(sizeof(ekind->tcstat[i].vscale_nhc),
1330 &(ekind->tcstat[i].vscale_nhc), cr);
1332 gmx_bcast(DIM*DIM*sizeof(ekind->ekin[0][0]),
1333 ekind->ekin[0], cr);
1335 gmx_bcast(sizeof(ekind->dekindl), &ekind->dekindl, cr);
1336 gmx_bcast(sizeof(ekind->cosacc.mvcos), &ekind->cosacc.mvcos, cr);
1340 void update_tcouple(int64_t step,
1341 const t_inputrec *inputrec,
1343 gmx_ekindata_t *ekind,
1344 const t_extmass *MassQ,
1345 const t_mdatoms *md)
1348 bool doTemperatureCoupling = false;
1350 /* if using vv with trotter decomposition methods, we do this elsewhere in the code */
1351 if (!(EI_VV(inputrec->eI) &&
1352 (inputrecNvtTrotter(inputrec) || inputrecNptTrotter(inputrec) || inputrecNphTrotter(inputrec))))
1354 doTemperatureCoupling = isTemperatureCouplingStep(step, inputrec);
1357 if (doTemperatureCoupling)
1359 real dttc = inputrec->nsttcouple*inputrec->delta_t;
1361 switch (inputrec->etc)
1366 berendsen_tcoupl(inputrec, ekind, dttc, state->therm_integral);
1369 nosehoover_tcoupl(&(inputrec->opts), ekind, dttc,
1370 state->nosehoover_xi.data(), state->nosehoover_vxi.data(), MassQ);
1373 vrescale_tcoupl(inputrec, step, ekind, dttc,
1374 state->therm_integral.data());
1377 /* rescale in place here */
1378 if (EI_VV(inputrec->eI))
1380 rescale_velocities(ekind, md, 0, md->homenr, state->v.rvec_array());
1385 /* Set the T scaling lambda to 1 to have no scaling */
1386 for (int i = 0; (i < inputrec->opts.ngtc); i++)
1388 ekind->tcstat[i].lambda = 1.0;
1393 /*! \brief Computes the atom range for a thread to operate on, ensured SIMD aligned ranges
1395 * \param[in] numThreads The number of threads to divide atoms over
1396 * \param[in] threadIndex The thread to get the range for
1397 * \param[in] numAtoms The total number of atoms (on this rank)
1398 * \param[out] startAtom The start of the atom range
1399 * \param[out] endAtom The end of the atom range, note that this is in general not a multiple of the SIMD width
1401 static void getThreadAtomRange(int numThreads, int threadIndex, int numAtoms,
1402 int *startAtom, int *endAtom)
1404 #if GMX_HAVE_SIMD_UPDATE
1405 constexpr int blockSize = GMX_SIMD_REAL_WIDTH;
1407 constexpr int blockSize = 1;
1409 int numBlocks = (numAtoms + blockSize - 1)/blockSize;
1411 *startAtom = ((numBlocks* threadIndex )/numThreads)*blockSize;
1412 *endAtom = ((numBlocks*(threadIndex + 1))/numThreads)*blockSize;
1413 if (threadIndex == numThreads - 1)
1415 *endAtom = numAtoms;
1419 void update_pcouple_before_coordinates(FILE *fplog,
1421 const t_inputrec *inputrec,
1423 matrix parrinellorahmanMu,
1427 /* Berendsen P-coupling is completely handled after the coordinate update.
1428 * Trotter P-coupling is handled by separate calls to trotter_update().
1430 if (inputrec->epc == epcPARRINELLORAHMAN &&
1431 isPressureCouplingStep(step, inputrec))
1433 real dtpc = inputrec->nstpcouple*inputrec->delta_t;
1435 parrinellorahman_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev,
1436 state->box, state->box_rel, state->boxv,
1437 M, parrinellorahmanMu, bInitStep);
1441 void constrain_velocities(int64_t step,
1442 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1445 gmx::Constraints *constr,
1457 * APPLY CONSTRAINTS:
1464 /* clear out constraints before applying */
1465 clear_mat(vir_part);
1467 /* Constrain the coordinates upd->xp */
1468 constr->apply(do_log, do_ene,
1470 state->x.rvec_array(), state->v.rvec_array(), state->v.rvec_array(),
1472 state->lambda[efptBONDED], dvdlambda,
1473 nullptr, bCalcVir ? &vir_con : nullptr, ConstraintVariable::Velocities);
1477 m_add(vir_part, vir_con, vir_part);
1482 void constrain_coordinates(int64_t step,
1483 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1487 gmx::Constraints *constr,
1500 /* clear out constraints before applying */
1501 clear_mat(vir_part);
1503 /* Constrain the coordinates upd->xp */
1504 constr->apply(do_log, do_ene,
1506 state->x.rvec_array(), upd->xp.rvec_array(), nullptr,
1508 state->lambda[efptBONDED], dvdlambda,
1509 as_rvec_array(state->v.data()), bCalcVir ? &vir_con : nullptr, ConstraintVariable::Positions);
1513 m_add(vir_part, vir_con, vir_part);
1519 update_sd_second_half(int64_t step,
1520 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1521 const t_inputrec *inputrec, /* input record and box stuff */
1522 const t_mdatoms *md,
1524 const t_commrec *cr,
1526 gmx_wallcycle_t wcycle,
1528 gmx::Constraints *constr,
1536 if (inputrec->eI == eiSD1)
1539 int homenr = md->homenr;
1541 /* Cast delta_t from double to real to make the integrators faster.
1542 * The only reason for having delta_t double is to get accurate values
1543 * for t=delta_t*step when step is larger than float precision.
1544 * For integration dt the accuracy of real suffices, since with
1545 * integral += dt*integrand the increment is nearly always (much) smaller
1546 * than the integral (and the integrand has real precision).
1548 real dt = inputrec->delta_t;
1550 wallcycle_start(wcycle, ewcUPDATE);
1552 nth = gmx_omp_nthreads_get(emntUpdate);
1554 #pragma omp parallel for num_threads(nth) schedule(static)
1555 for (th = 0; th < nth; th++)
1559 int start_th, end_th;
1560 getThreadAtomRange(nth, th, homenr, &start_th, &end_th);
1562 doSDUpdateGeneral<SDUpdate::FrictionAndNoiseOnly>
1564 start_th, end_th, dt,
1565 inputrec->opts.acc, inputrec->opts.nFreeze,
1566 md->invmass, md->ptype,
1567 md->cFREEZE, nullptr, md->cTC,
1568 state->x.rvec_array(), upd->xp.rvec_array(),
1569 state->v.rvec_array(), nullptr,
1570 step, inputrec->ld_seed,
1571 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1573 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1575 inc_nrnb(nrnb, eNR_UPDATE, homenr);
1576 wallcycle_stop(wcycle, ewcUPDATE);
1578 /* Constrain the coordinates upd->xp for half a time step */
1579 constr->apply(do_log, do_ene,
1581 state->x.rvec_array(), upd->xp.rvec_array(), nullptr,
1583 state->lambda[efptBONDED], dvdlambda,
1584 as_rvec_array(state->v.data()), nullptr, ConstraintVariable::Positions);
1588 void finish_update(const t_inputrec *inputrec, /* input record and box stuff */
1589 const t_mdatoms *md,
1591 const t_graph *graph,
1593 gmx_wallcycle_t wcycle,
1595 const gmx::Constraints *constr)
1597 int homenr = md->homenr;
1599 /* We must always unshift after updating coordinates; if we did not shake
1600 x was shifted in do_force */
1602 /* NOTE Currently we always integrate to a temporary buffer and
1603 * then copy the results back. */
1605 wallcycle_start_nocount(wcycle, ewcUPDATE);
1607 if (md->cFREEZE != nullptr && constr != nullptr)
1609 /* If we have atoms that are frozen along some, but not all
1610 * dimensions, then any constraints will have moved them also along
1611 * the frozen dimensions. To freeze such degrees of freedom
1612 * we copy them back here to later copy them forward. It would
1613 * be more elegant and slightly more efficient to copies zero
1614 * times instead of twice, but the graph case below prevents this.
1616 const ivec *nFreeze = inputrec->opts.nFreeze;
1617 bool partialFreezeAndConstraints = false;
1618 for (int g = 0; g < inputrec->opts.ngfrz; g++)
1620 int numFreezeDim = nFreeze[g][XX] + nFreeze[g][YY] + nFreeze[g][ZZ];
1621 if (numFreezeDim > 0 && numFreezeDim < 3)
1623 partialFreezeAndConstraints = true;
1626 if (partialFreezeAndConstraints)
1628 auto xp = makeArrayRef(upd->xp).subArray(0, homenr);
1629 auto x = makeConstArrayRef(state->x).subArray(0, homenr);
1630 for (int i = 0; i < homenr; i++)
1632 int g = md->cFREEZE[i];
1634 for (int d = 0; d < DIM; d++)
1645 if (graph && (graph->nnodes > 0))
1647 unshift_x(graph, state->box, state->x.rvec_array(), upd->xp.rvec_array());
1648 if (TRICLINIC(state->box))
1650 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
1654 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
1659 auto xp = makeConstArrayRef(upd->xp).subArray(0, homenr);
1660 auto x = makeArrayRef(state->x).subArray(0, homenr);
1661 #ifndef __clang_analyzer__
1662 int gmx_unused nth = gmx_omp_nthreads_get(emntUpdate);
1664 #pragma omp parallel for num_threads(nth) schedule(static)
1665 for (int i = 0; i < homenr; i++)
1667 // Trivial statement, does not throw
1671 wallcycle_stop(wcycle, ewcUPDATE);
1673 /* ############# END the update of velocities and positions ######### */
1676 void update_pcouple_after_coordinates(FILE *fplog,
1678 const t_inputrec *inputrec,
1679 const t_mdatoms *md,
1680 const matrix pressure,
1681 const matrix forceVirial,
1682 const matrix constraintVirial,
1683 const matrix parrinellorahmanMu,
1689 int homenr = md->homenr;
1691 /* Cast to real for faster code, no loss in precision (see comment above) */
1692 real dt = inputrec->delta_t;
1695 /* now update boxes */
1696 switch (inputrec->epc)
1700 case (epcBERENDSEN):
1701 if (isPressureCouplingStep(step, inputrec))
1703 real dtpc = inputrec->nstpcouple*dt;
1705 berendsen_pcoupl(fplog, step, inputrec, dtpc,
1706 pressure, state->box,
1707 forceVirial, constraintVirial,
1708 mu, &state->baros_integral);
1709 berendsen_pscale(inputrec, mu, state->box, state->box_rel,
1710 start, homenr, state->x.rvec_array(),
1714 case (epcPARRINELLORAHMAN):
1715 if (isPressureCouplingStep(step, inputrec))
1717 /* The box velocities were updated in do_pr_pcoupl,
1718 * but we dont change the box vectors until we get here
1719 * since we need to be able to shift/unshift above.
1721 real dtpc = inputrec->nstpcouple*dt;
1722 for (int i = 0; i < DIM; i++)
1724 for (int m = 0; m <= i; m++)
1726 state->box[i][m] += dtpc*state->boxv[i][m];
1729 preserve_box_shape(inputrec, state->box_rel, state->box);
1731 /* Scale the coordinates */
1732 auto x = state->x.rvec_array();
1733 for (int n = start; n < start + homenr; n++)
1735 tmvmul_ur0(parrinellorahmanMu, x[n], x[n]);
1740 switch (inputrec->epct)
1742 case (epctISOTROPIC):
1743 /* DIM * eta = ln V. so DIM*eta_new = DIM*eta_old + DIM*dt*veta =>
1744 ln V_new = ln V_old + 3*dt*veta => V_new = V_old*exp(3*dt*veta) =>
1745 Side length scales as exp(veta*dt) */
1747 msmul(state->box, std::exp(state->veta*dt), state->box);
1749 /* Relate veta to boxv. veta = d(eta)/dT = (1/DIM)*1/V dV/dT.
1750 o If we assume isotropic scaling, and box length scaling
1751 factor L, then V = L^DIM (det(M)). So dV/dt = DIM
1752 L^(DIM-1) dL/dt det(M), and veta = (1/L) dL/dt. The
1753 determinant of B is L^DIM det(M), and the determinant
1754 of dB/dt is (dL/dT)^DIM det (M). veta will be
1755 (det(dB/dT)/det(B))^(1/3). Then since M =
1756 B_new*(vol_new)^(1/3), dB/dT_new = (veta_new)*B(new). */
1758 msmul(state->box, state->veta, state->boxv);
1770 auto localX = makeArrayRef(state->x).subArray(start, homenr);
1771 upd->deform->apply(localX, state->box, step);
1775 void update_coords(int64_t step,
1776 const t_inputrec *inputrec, /* input record and box stuff */
1777 const t_mdatoms *md,
1779 gmx::ArrayRefWithPadding<gmx::RVec> f,
1780 const t_fcdata *fcd,
1781 const gmx_ekindata_t *ekind,
1785 const t_commrec *cr, /* these shouldn't be here -- need to think about it */
1786 const gmx::Constraints *constr)
1788 gmx_bool bDoConstr = (nullptr != constr);
1790 /* Running the velocity half does nothing except for velocity verlet */
1791 if ((UpdatePart == etrtVELOCITY1 || UpdatePart == etrtVELOCITY2) &&
1792 !EI_VV(inputrec->eI))
1794 gmx_incons("update_coords called for velocity without VV integrator");
1797 int homenr = md->homenr;
1799 /* Cast to real for faster code, no loss in precision (see comment above) */
1800 real dt = inputrec->delta_t;
1802 /* We need to update the NMR restraint history when time averaging is used */
1803 if (state->flags & (1<<estDISRE_RM3TAV))
1805 update_disres_history(fcd, &state->hist);
1807 if (state->flags & (1<<estORIRE_DTAV))
1809 update_orires_history(fcd, &state->hist);
1812 /* ############# START The update of velocities and positions ######### */
1813 int nth = gmx_omp_nthreads_get(emntUpdate);
1815 #pragma omp parallel for num_threads(nth) schedule(static)
1816 for (int th = 0; th < nth; th++)
1820 int start_th, end_th;
1821 getThreadAtomRange(nth, th, homenr, &start_th, &end_th);
1823 const rvec *x_rvec = state->x.rvec_array();
1824 rvec *xp_rvec = upd->xp.rvec_array();
1825 rvec *v_rvec = state->v.rvec_array();
1826 const rvec *f_rvec = as_rvec_array(f.unpaddedArrayRef().data());
1828 switch (inputrec->eI)
1831 do_update_md(start_th, end_th, step, dt,
1832 inputrec, md, ekind, state->box,
1833 x_rvec, xp_rvec, v_rvec, f_rvec,
1834 state->nosehoover_vxi.data(), M);
1839 // With constraints, the SD update is done in 2 parts
1840 doSDUpdateGeneral<SDUpdate::ForcesOnly>
1842 start_th, end_th, dt,
1843 inputrec->opts.acc, inputrec->opts.nFreeze,
1844 md->invmass, md->ptype,
1845 md->cFREEZE, md->cACC, nullptr,
1846 x_rvec, xp_rvec, v_rvec, f_rvec,
1847 step, inputrec->ld_seed, nullptr);
1851 doSDUpdateGeneral<SDUpdate::Combined>
1853 start_th, end_th, dt,
1854 inputrec->opts.acc, inputrec->opts.nFreeze,
1855 md->invmass, md->ptype,
1856 md->cFREEZE, md->cACC, md->cTC,
1857 x_rvec, xp_rvec, v_rvec, f_rvec,
1858 step, inputrec->ld_seed,
1859 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1863 do_update_bd(start_th, end_th, dt,
1864 inputrec->opts.nFreeze, md->invmass, md->ptype,
1865 md->cFREEZE, md->cTC,
1866 x_rvec, xp_rvec, v_rvec, f_rvec,
1868 upd->sd->bd_rf.data(),
1869 step, inputrec->ld_seed, DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1874 gmx_bool bExtended = (inputrec->etc == etcNOSEHOOVER ||
1875 inputrec->epc == epcPARRINELLORAHMAN ||
1876 inputrec->epc == epcMTTK);
1878 /* assuming barostat coupled to group 0 */
1879 real alpha = 1.0 + DIM/static_cast<real>(inputrec->opts.nrdf[0]);
1884 do_update_vv_vel(start_th, end_th, dt,
1885 inputrec->opts.acc, inputrec->opts.nFreeze,
1886 md->invmass, md->ptype,
1887 md->cFREEZE, md->cACC,
1889 bExtended, state->veta, alpha);
1892 do_update_vv_pos(start_th, end_th, dt,
1893 inputrec->opts.nFreeze,
1894 md->ptype, md->cFREEZE,
1895 x_rvec, xp_rvec, v_rvec,
1896 bExtended, state->veta);
1902 gmx_fatal(FARGS, "Don't know how to update coordinates");
1905 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1910 extern gmx_bool update_randomize_velocities(const t_inputrec *ir, int64_t step, const t_commrec *cr,
1911 const t_mdatoms *md,
1912 gmx::ArrayRef<gmx::RVec> v,
1913 const gmx_update_t *upd,
1914 const gmx::Constraints *constr)
1917 real rate = (ir->delta_t)/ir->opts.tau_t[0];
1919 if (ir->etc == etcANDERSEN && constr != nullptr)
1921 /* Currently, Andersen thermostat does not support constrained
1922 systems. Functionality exists in the andersen_tcoupl
1923 function in GROMACS 4.5.7 to allow this combination. That
1924 code could be ported to the current random-number
1925 generation approach, but has not yet been done because of
1926 lack of time and resources. */
1927 gmx_fatal(FARGS, "Normal Andersen is currently not supported with constraints, use massive Andersen instead");
1930 /* proceed with andersen if 1) it's fixed probability per
1931 particle andersen or 2) it's massive andersen and it's tau_t/dt */
1932 if ((ir->etc == etcANDERSEN) || do_per_step(step, roundToInt(1.0/rate)))
1934 andersen_tcoupl(ir, step, cr, md, v, rate,
1935 upd->sd->randomize_group,