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, 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/domdec/domdec_struct.h"
47 #include "gromacs/fileio/confio.h"
48 #include "gromacs/gmxlib/network.h"
49 #include "gromacs/gmxlib/nrnb.h"
50 #include "gromacs/listed-forces/disre.h"
51 #include "gromacs/listed-forces/orires.h"
52 #include "gromacs/math/functions.h"
53 #include "gromacs/math/invertmatrix.h"
54 #include "gromacs/math/paddedvector.h"
55 #include "gromacs/math/units.h"
56 #include "gromacs/math/vec.h"
57 #include "gromacs/math/vecdump.h"
58 #include "gromacs/mdlib/boxdeformation.h"
59 #include "gromacs/mdlib/constr.h"
60 #include "gromacs/mdlib/gmx_omp_nthreads.h"
61 #include "gromacs/mdlib/mdrun.h"
62 #include "gromacs/mdlib/sim_util.h"
63 #include "gromacs/mdlib/tgroup.h"
64 #include "gromacs/mdtypes/commrec.h"
65 #include "gromacs/mdtypes/group.h"
66 #include "gromacs/mdtypes/inputrec.h"
67 #include "gromacs/mdtypes/md_enums.h"
68 #include "gromacs/mdtypes/state.h"
69 #include "gromacs/pbcutil/boxutilities.h"
70 #include "gromacs/pbcutil/mshift.h"
71 #include "gromacs/pbcutil/pbc.h"
72 #include "gromacs/pulling/pull.h"
73 #include "gromacs/random/tabulatednormaldistribution.h"
74 #include "gromacs/random/threefry.h"
75 #include "gromacs/simd/simd.h"
76 #include "gromacs/timing/wallcycle.h"
77 #include "gromacs/topology/atoms.h"
78 #include "gromacs/utility/exceptions.h"
79 #include "gromacs/utility/fatalerror.h"
80 #include "gromacs/utility/futil.h"
81 #include "gromacs/utility/gmxassert.h"
82 #include "gromacs/utility/gmxomp.h"
83 #include "gromacs/utility/smalloc.h"
85 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
100 gmx_sd_sigma_t *sdsig;
101 /* andersen temperature control stuff */
102 gmx_bool *randomize_group;
109 /* xprime for constraint algorithms */
110 PaddedVector<gmx::RVec> xp;
112 /* Variables for the deform algorithm */
113 int64_t deformref_step;
114 matrix deformref_box;
116 //! Box deformation handler (or nullptr if inactive).
117 gmx::BoxDeformation *deform;
120 static bool isTemperatureCouplingStep(int64_t step, const t_inputrec *ir)
122 /* We should only couple after a step where energies were determined (for leapfrog versions)
123 or the step energies are determined, for velocity verlet versions */
133 return ir->etc != etcNO &&
134 (ir->nsttcouple == 1 ||
135 do_per_step(step + ir->nsttcouple - offset, ir->nsttcouple));
138 static bool isPressureCouplingStep(int64_t step, const t_inputrec *ir)
140 GMX_ASSERT(ir->epc != epcMTTK, "MTTK pressure coupling is not handled here");
143 if (ir->epc == epcBERENDSEN)
151 /* We should only couple after a step where pressures were determined */
152 return ir->epc != etcNO &&
153 (ir->nstpcouple == 1 ||
154 do_per_step(step + ir->nstpcouple - offset, ir->nstpcouple));
157 /*! \brief Sets the velocities of virtual sites to zero */
158 static void clearVsiteVelocities(int start,
160 const unsigned short *particleType,
161 rvec * gmx_restrict v)
163 for (int a = start; a < nrend; a++)
165 if (particleType[a] == eptVSite)
172 /*! \brief Sets the number of different temperature coupling values */
173 enum class NumTempScaleValues
175 single, //!< Single T-scaling value (either one group or all values =1)
176 multiple //!< Multiple T-scaling values, need to use T-group indices
179 /*! \brief Sets if to apply no or diagonal Parrinello-Rahman pressure scaling
181 * Note that this enum is only used in updateMDLeapfrogSimple(), which does
182 * not handle fully anistropic Parrinello-Rahman scaling, so we only have
183 * options \p no and \p diagonal here and no anistropic option.
185 enum class ApplyParrinelloRahmanVScaling
187 no, //!< Do not apply velocity scaling (not a PR-coupling run or step)
188 diagonal //!< Apply velocity scaling using a diagonal matrix
191 /*! \brief Integrate using leap-frog with T-scaling and optionally diagonal Parrinello-Rahman p-coupling
193 * \tparam numTempScaleValues The number of different T-couple values
194 * \tparam applyPRVScaling Apply Parrinello-Rahman velocity scaling
195 * \param[in] start Index of first atom to update
196 * \param[in] nrend Last atom to update: \p nrend - 1
197 * \param[in] dt The time step
198 * \param[in] dtPressureCouple Time step for pressure coupling
199 * \param[in] invMassPerDim 1/mass per atom and dimension
200 * \param[in] tcstat Temperature coupling information
201 * \param[in] cTC T-coupling group index per atom
202 * \param[in] pRVScaleMatrixDiagonal Parrinello-Rahman v-scale matrix diagonal
203 * \param[in] x Input coordinates
204 * \param[out] xprime Updated coordinates
205 * \param[inout] v Velocities
206 * \param[in] f Forces
208 * We expect this template to get good SIMD acceleration by most compilers,
209 * unlike the more complex general template.
210 * Note that we might get even better SIMD acceleration when we introduce
211 * aligned (and padded) memory, possibly with some hints for the compilers.
213 template<NumTempScaleValues numTempScaleValues,
214 ApplyParrinelloRahmanVScaling applyPRVScaling>
216 updateMDLeapfrogSimple(int start,
219 real dtPressureCouple,
220 const rvec * gmx_restrict invMassPerDim,
221 const t_grp_tcstat * tcstat,
222 const unsigned short * cTC,
223 const rvec pRVScaleMatrixDiagonal,
224 const rvec * gmx_restrict x,
225 rvec * gmx_restrict xprime,
226 rvec * gmx_restrict v,
227 const rvec * gmx_restrict f)
231 if (numTempScaleValues == NumTempScaleValues::single)
233 lambdaGroup = tcstat[0].lambda;
236 for (int a = start; a < nrend; a++)
238 if (numTempScaleValues == NumTempScaleValues::multiple)
240 lambdaGroup = tcstat[cTC[a]].lambda;
243 for (int d = 0; d < DIM; d++)
245 /* Note that using rvec invMassPerDim results in more efficient
246 * SIMD code, but this increases the cache pressure.
247 * For large systems with PME on the CPU this slows down the
248 * (then already slow) update by 20%. If all data remains in cache,
249 * using rvec is much faster.
251 real vNew = lambdaGroup*v[a][d] + f[a][d]*invMassPerDim[a][d]*dt;
253 if (applyPRVScaling == ApplyParrinelloRahmanVScaling::diagonal)
255 vNew -= dtPressureCouple*pRVScaleMatrixDiagonal[d]*v[a][d];
258 xprime[a][d] = x[a][d] + vNew*dt;
263 #if GMX_SIMD && GMX_SIMD_HAVE_REAL
264 #define GMX_HAVE_SIMD_UPDATE 1
266 #define GMX_HAVE_SIMD_UPDATE 0
269 #if GMX_HAVE_SIMD_UPDATE
271 /*! \brief Load (aligned) the contents of GMX_SIMD_REAL_WIDTH rvec elements sequentially into 3 SIMD registers
273 * The loaded output is:
274 * \p r0: { r[index][XX], r[index][YY], ... }
276 * \p r2: { ..., r[index+GMX_SIMD_REAL_WIDTH-1][YY], r[index+GMX_SIMD_REAL_WIDTH-1][ZZ] }
278 * \param[in] r Real to an rvec array, has to be aligned to SIMD register width
279 * \param[in] index Index of the first rvec triplet of reals to load
280 * \param[out] r0 Pointer to first SIMD register
281 * \param[out] r1 Pointer to second SIMD register
282 * \param[out] r2 Pointer to third SIMD register
284 static inline void simdLoadRvecs(const rvec *r,
290 const real *realPtr = r[index];
292 GMX_ASSERT(isSimdAligned(realPtr), "Pointer should be SIMD aligned");
294 *r0 = simdLoad(realPtr + 0*GMX_SIMD_REAL_WIDTH);
295 *r1 = simdLoad(realPtr + 1*GMX_SIMD_REAL_WIDTH);
296 *r2 = simdLoad(realPtr + 2*GMX_SIMD_REAL_WIDTH);
299 /*! \brief Store (aligned) 3 SIMD registers sequentially to GMX_SIMD_REAL_WIDTH rvec elements
301 * The stored output is:
302 * \p r[index] = { { r0[0], r0[1], ... }
304 * \p r[index+GMX_SIMD_REAL_WIDTH-1] = { ... , r2[GMX_SIMD_REAL_WIDTH-2], r2[GMX_SIMD_REAL_WIDTH-1] }
306 * \param[out] r Pointer 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 store to
308 * \param[in] r0 First SIMD register
309 * \param[in] r1 Second SIMD register
310 * \param[in] r2 Third SIMD register
312 static inline void simdStoreRvecs(rvec *r,
318 real *realPtr = r[index];
320 GMX_ASSERT(isSimdAligned(realPtr), "Pointer should be SIMD aligned");
322 store(realPtr + 0*GMX_SIMD_REAL_WIDTH, r0);
323 store(realPtr + 1*GMX_SIMD_REAL_WIDTH, r1);
324 store(realPtr + 2*GMX_SIMD_REAL_WIDTH, r2);
327 /*! \brief Integrate using leap-frog with single group T-scaling and SIMD
329 * \param[in] start Index of first atom to update
330 * \param[in] nrend Last atom to update: \p nrend - 1
331 * \param[in] dt The time step
332 * \param[in] invMass 1/mass per atom
333 * \param[in] tcstat Temperature coupling information
334 * \param[in] x Input coordinates
335 * \param[out] xprime Updated coordinates
336 * \param[inout] v Velocities
337 * \param[in] f Forces
340 updateMDLeapfrogSimpleSimd(int start,
343 const real * gmx_restrict invMass,
344 const t_grp_tcstat * tcstat,
345 const rvec * gmx_restrict x,
346 rvec * gmx_restrict xprime,
347 rvec * gmx_restrict v,
348 const rvec * gmx_restrict f)
350 SimdReal timestep(dt);
351 SimdReal lambdaSystem(tcstat[0].lambda);
353 /* We declare variables here, since code is often slower when declaring them inside the loop */
355 /* Note: We should implement a proper PaddedVector, so we don't need this check */
356 GMX_ASSERT(isSimdAligned(invMass), "invMass should be aligned");
358 for (int a = start; a < nrend; a += GMX_SIMD_REAL_WIDTH)
360 SimdReal invMass0, invMass1, invMass2;
361 expandScalarsToTriplets(simdLoad(invMass + a),
362 &invMass0, &invMass1, &invMass2);
366 simdLoadRvecs(v, a, &v0, &v1, &v2);
367 simdLoadRvecs(f, a, &f0, &f1, &f2);
369 v0 = fma(f0*invMass0, timestep, lambdaSystem*v0);
370 v1 = fma(f1*invMass1, timestep, lambdaSystem*v1);
371 v2 = fma(f2*invMass2, timestep, lambdaSystem*v2);
373 simdStoreRvecs(v, a, v0, v1, v2);
376 simdLoadRvecs(x, a, &x0, &x1, &x2);
378 SimdReal xprime0 = fma(v0, timestep, x0);
379 SimdReal xprime1 = fma(v1, timestep, x1);
380 SimdReal xprime2 = fma(v2, timestep, x2);
382 simdStoreRvecs(xprime, a, xprime0, xprime1, xprime2);
386 #endif // GMX_HAVE_SIMD_UPDATE
388 /*! \brief Sets the NEMD acceleration type */
389 enum class AccelerationType
394 /*! \brief Integrate using leap-frog with support for everything
396 * \tparam accelerationType Type of NEMD acceleration
397 * \param[in] start Index of first atom to update
398 * \param[in] nrend Last atom to update: \p nrend - 1
399 * \param[in] doNoseHoover If to apply Nose-Hoover T-coupling
400 * \param[in] dt The time step
401 * \param[in] dtPressureCouple Time step for pressure coupling, is 0 when no pressure coupling should be applied at this step
402 * \param[in] ir The input parameter record
403 * \param[in] md Atom properties
404 * \param[in] ekind Kinetic energy data
405 * \param[in] box The box dimensions
406 * \param[in] x Input coordinates
407 * \param[out] xprime Updated coordinates
408 * \param[inout] v Velocities
409 * \param[in] f Forces
410 * \param[in] nh_vxi Nose-Hoover velocity scaling factors
411 * \param[in] M Parrinello-Rahman scaling matrix
413 template<AccelerationType accelerationType>
415 updateMDLeapfrogGeneral(int start,
419 real dtPressureCouple,
420 const t_inputrec * ir,
421 const t_mdatoms * md,
422 const gmx_ekindata_t * ekind,
424 const rvec * gmx_restrict x,
425 rvec * gmx_restrict xprime,
426 rvec * gmx_restrict v,
427 const rvec * gmx_restrict f,
428 const double * gmx_restrict nh_vxi,
431 /* This is a version of the leap-frog integrator that supports
432 * all combinations of T-coupling, P-coupling and NEMD.
433 * Nose-Hoover uses the reversible leap-frog integrator from
434 * Holian et al. Phys Rev E 52(3) : 2338, 1995
437 const unsigned short * cTC = md->cTC;
438 const t_grp_tcstat * tcstat = ekind->tcstat;
440 const unsigned short * cACC = md->cACC;
441 const rvec * accel = ir->opts.acc;
442 const t_grp_acc * grpstat = ekind->grpstat;
444 const rvec * gmx_restrict invMassPerDim = md->invMassPerDim;
446 /* Initialize group values, changed later when multiple groups are used */
451 for (int n = start; n < nrend; n++)
457 real lg = tcstat[gt].lambda;
460 real cosineZ, vCosine;
461 #ifdef __INTEL_COMPILER
462 #pragma warning( disable : 280 )
464 switch (accelerationType)
466 case AccelerationType::none:
467 copy_rvec(v[n], vRel);
469 case AccelerationType::group:
474 /* Avoid scaling the group velocity */
475 rvec_sub(v[n], grpstat[ga].u, vRel);
477 case AccelerationType::cosine:
478 cosineZ = std::cos(x[n][ZZ]*static_cast<real>(M_PI)/box[ZZ][ZZ]);
479 vCosine = cosineZ*ekind->cosacc.vcos;
480 /* Avoid scaling the cosine profile velocity */
481 copy_rvec(v[n], vRel);
488 /* Here we account for multiple time stepping, by increasing
489 * the Nose-Hoover correction by nsttcouple
491 factorNH = 0.5*ir->nsttcouple*dt*nh_vxi[gt];
494 for (int d = 0; d < DIM; d++)
497 (lg*vRel[d] + (f[n][d]*invMassPerDim[n][d]*dt - factorNH*vRel[d]
498 - dtPressureCouple*iprod(M[d], vRel)))/(1 + factorNH);
499 switch (accelerationType)
501 case AccelerationType::none:
503 case AccelerationType::group:
504 /* Add back the mean velocity and apply acceleration */
505 vNew += grpstat[ga].u[d] + accel[ga][d]*dt;
507 case AccelerationType::cosine:
510 /* Add back the mean velocity and apply acceleration */
511 vNew += vCosine + cosineZ*ekind->cosacc.cos_accel*dt;
516 xprime[n][d] = x[n][d] + vNew*dt;
521 /*! \brief Handles the Leap-frog MD x and v integration */
522 static void do_update_md(int start,
526 const t_inputrec * ir,
527 const t_mdatoms * md,
528 const gmx_ekindata_t * ekind,
530 const rvec * gmx_restrict x,
531 rvec * gmx_restrict xprime,
532 rvec * gmx_restrict v,
533 const rvec * gmx_restrict f,
534 const double * gmx_restrict nh_vxi,
537 GMX_ASSERT(nrend == start || xprime != x, "For SIMD optimization certain compilers need to have xprime != x");
539 /* Note: Berendsen pressure scaling is handled after do_update_md() */
540 bool doTempCouple = isTemperatureCouplingStep(step, ir);
541 bool doNoseHoover = (ir->etc == etcNOSEHOOVER && doTempCouple);
542 bool doParrinelloRahman = (ir->epc == epcPARRINELLORAHMAN && isPressureCouplingStep(step, ir));
543 bool doPROffDiagonal = (doParrinelloRahman && (M[YY][XX] != 0 || M[ZZ][XX] != 0 || M[ZZ][YY] != 0));
545 real dtPressureCouple = (doParrinelloRahman ? ir->nstpcouple*dt : 0);
547 /* NEMD (also cosine) acceleration is applied in updateMDLeapFrogGeneral */
548 bool doAcceleration = (ekind->bNEMD || ekind->cosacc.cos_accel != 0);
550 if (doNoseHoover || doPROffDiagonal || doAcceleration)
553 if (!doParrinelloRahman)
555 /* We should not apply PR scaling at this step */
565 updateMDLeapfrogGeneral<AccelerationType::none>
566 (start, nrend, doNoseHoover, dt, dtPressureCouple,
567 ir, md, ekind, box, x, xprime, v, f, nh_vxi, stepM);
569 else if (ekind->bNEMD)
571 updateMDLeapfrogGeneral<AccelerationType::group>
572 (start, nrend, doNoseHoover, dt, dtPressureCouple,
573 ir, md, ekind, box, x, xprime, v, f, nh_vxi, stepM);
577 updateMDLeapfrogGeneral<AccelerationType::cosine>
578 (start, nrend, doNoseHoover, dt, dtPressureCouple,
579 ir, md, ekind, box, x, xprime, v, f, nh_vxi, stepM);
584 /* Use a simple and thus more efficient integration loop. */
585 /* The simple loop does not check for particle type (so it can
586 * be vectorized), which means we need to clear the velocities
587 * of virtual sites in advance, when present. Note that vsite
588 * velocities are computed after update and constraints from
589 * their displacement.
593 /* Note: The overhead of this loop is completely neligible */
594 clearVsiteVelocities(start, nrend, md->ptype, v);
597 /* We determine if we have a single T-coupling lambda value for all
598 * atoms. That allows for better SIMD acceleration in the template.
599 * If we do not do temperature coupling (in the run or this step),
600 * all scaling values are 1, so we effectively have a single value.
602 bool haveSingleTempScaleValue = (!doTempCouple || ekind->ngtc == 1);
604 /* Extract some pointers needed by all cases */
605 const unsigned short *cTC = md->cTC;
606 const t_grp_tcstat *tcstat = ekind->tcstat;
607 const rvec *invMassPerDim = md->invMassPerDim;
609 if (doParrinelloRahman)
611 GMX_ASSERT(!doPROffDiagonal, "updateMDLeapfrogSimple only support diagonal Parrinello-Rahman scaling matrices");
614 for (int d = 0; d < DIM; d++)
619 if (haveSingleTempScaleValue)
621 updateMDLeapfrogSimple
622 <NumTempScaleValues::single,
623 ApplyParrinelloRahmanVScaling::diagonal>
624 (start, nrend, dt, dtPressureCouple,
625 invMassPerDim, tcstat, cTC, diagM, x, xprime, v, f);
629 updateMDLeapfrogSimple
630 <NumTempScaleValues::multiple,
631 ApplyParrinelloRahmanVScaling::diagonal>
632 (start, nrend, dt, dtPressureCouple,
633 invMassPerDim, tcstat, cTC, diagM, x, xprime, v, f);
638 if (haveSingleTempScaleValue)
640 /* Note that modern compilers are pretty good at vectorizing
641 * updateMDLeapfrogSimple(). But the SIMD version will still
642 * be faster because invMass lowers the cache pressure
643 * compared to invMassPerDim.
645 #if GMX_HAVE_SIMD_UPDATE
646 /* Check if we can use invmass instead of invMassPerDim */
647 if (!md->havePartiallyFrozenAtoms)
649 updateMDLeapfrogSimpleSimd
650 (start, nrend, dt, md->invmass, tcstat, x, xprime, v, f);
655 updateMDLeapfrogSimple
656 <NumTempScaleValues::single,
657 ApplyParrinelloRahmanVScaling::no>
658 (start, nrend, dt, dtPressureCouple,
659 invMassPerDim, tcstat, cTC, nullptr, x, xprime, v, f);
664 updateMDLeapfrogSimple
665 <NumTempScaleValues::multiple,
666 ApplyParrinelloRahmanVScaling::no>
667 (start, nrend, dt, dtPressureCouple,
668 invMassPerDim, tcstat, cTC, nullptr, x, xprime, v, f);
674 static void do_update_vv_vel(int start, int nrend, real dt,
675 const rvec accel[], const ivec nFreeze[], const real invmass[],
676 const unsigned short ptype[], const unsigned short cFREEZE[],
677 const unsigned short cACC[], rvec v[], const rvec f[],
678 gmx_bool bExtended, real veta, real alpha)
686 g = 0.25*dt*veta*alpha;
688 mv2 = gmx::series_sinhx(g);
695 for (n = start; n < nrend; n++)
697 real w_dt = invmass[n]*dt;
707 for (d = 0; d < DIM; d++)
709 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
711 v[n][d] = mv1*(mv1*v[n][d] + 0.5*(w_dt*mv2*f[n][d]))+0.5*accel[ga][d]*dt;
719 } /* do_update_vv_vel */
721 static void do_update_vv_pos(int start, int nrend, real dt,
722 const ivec nFreeze[],
723 const unsigned short ptype[], const unsigned short cFREEZE[],
724 const rvec x[], rvec xprime[], const rvec v[],
725 gmx_bool bExtended, real veta)
731 /* Would it make more sense if Parrinello-Rahman was put here? */
736 mr2 = gmx::series_sinhx(g);
744 for (n = start; n < nrend; n++)
752 for (d = 0; d < DIM; d++)
754 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
756 xprime[n][d] = mr1*(mr1*x[n][d]+mr2*dt*v[n][d]);
760 xprime[n][d] = x[n][d];
764 } /* do_update_vv_pos */
766 static gmx_stochd_t *init_stochd(const t_inputrec *ir)
772 const t_grpopts *opts = &ir->opts;
773 int ngtc = opts->ngtc;
777 snew(sd->bd_rf, ngtc);
779 else if (EI_SD(ir->eI))
782 snew(sd->sdsig, ngtc);
784 gmx_sd_const_t *sdc = sd->sdc;
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 snew(sd->randomize_group, ngtc);
802 snew(sd->boltzfac, 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 sd->randomize_group[gt] = TRUE;
813 sd->boltzfac[gt] = BOLTZ*opts->ref_t[gt];
817 sd->randomize_group[gt] = FALSE;
825 void update_temperature_constants(gmx_update_t *upd, const t_inputrec *ir)
829 if (ir->bd_fric != 0)
831 for (int gt = 0; gt < ir->opts.ngtc; gt++)
833 upd->sd->bd_rf[gt] = std::sqrt(2.0*BOLTZ*ir->opts.ref_t[gt]/(ir->bd_fric*ir->delta_t));
838 for (int gt = 0; gt < ir->opts.ngtc; gt++)
840 upd->sd->bd_rf[gt] = std::sqrt(2.0*BOLTZ*ir->opts.ref_t[gt]);
846 for (int gt = 0; gt < ir->opts.ngtc; gt++)
848 real kT = BOLTZ*ir->opts.ref_t[gt];
849 /* The mass is accounted for later, since this differs per atom */
850 upd->sd->sdsig[gt].V = std::sqrt(kT*(1 - upd->sd->sdc[gt].em*upd->sd->sdc[gt].em));
855 gmx_update_t *init_update(const t_inputrec *ir,
856 gmx::BoxDeformation *deform)
858 gmx_update_t *upd = new(gmx_update_t);
860 if (ir->eI == eiBD || EI_SD(ir->eI) || ir->etc == etcVRESCALE || ETC_ANDERSEN(ir->etc))
862 upd->sd = init_stochd(ir);
865 update_temperature_constants(upd, ir);
867 upd->xp.resizeWithPadding(0);
869 upd->deform = deform;
874 void update_realloc(gmx_update_t *upd, int natoms)
876 GMX_ASSERT(upd, "upd must be allocated before its fields can be reallocated");
878 upd->xp.resizeWithPadding(natoms);
881 /*! \brief Sets the SD update type */
882 enum class SDUpdate : int
884 ForcesOnly, FrictionAndNoiseOnly, Combined
887 /*! \brief SD integrator update
889 * Two phases are required in the general case of a constrained
890 * update, the first phase from the contribution of forces, before
891 * applying constraints, and then a second phase applying the friction
892 * and noise, and then further constraining. For details, see
895 * Without constraints, the two phases can be combined, for
898 * Thus three instantiations of this templated function will be made,
899 * two with only one contribution, and one with both contributions. */
900 template <SDUpdate updateType>
902 doSDUpdateGeneral(gmx_stochd_t *sd,
903 int start, int nrend, real dt,
904 const rvec accel[], const ivec nFreeze[],
905 const real invmass[], const unsigned short ptype[],
906 const unsigned short cFREEZE[], const unsigned short cACC[],
907 const unsigned short cTC[],
908 const rvec x[], rvec xprime[], rvec v[], const rvec f[],
909 int64_t step, int seed, const int *gatindex)
911 // cTC, cACC and cFreeze can be nullptr any time, but various
912 // instantiations do not make sense with particular pointer
914 if (updateType == SDUpdate::ForcesOnly)
916 GMX_ASSERT(f != nullptr, "SD update with only forces requires forces");
917 GMX_ASSERT(cTC == nullptr, "SD update with only forces cannot handle temperature groups");
919 if (updateType == SDUpdate::FrictionAndNoiseOnly)
921 GMX_ASSERT(f == nullptr, "SD update with only noise cannot handle forces");
922 GMX_ASSERT(cACC == nullptr, "SD update with only noise cannot handle acceleration groups");
924 if (updateType == SDUpdate::Combined)
926 GMX_ASSERT(f != nullptr, "SD update with forces and noise requires forces");
929 // Even 0 bits internal counter gives 2x64 ints (more than enough for three table lookups)
930 gmx::ThreeFry2x64<0> rng(seed, gmx::RandomDomain::UpdateCoordinates);
931 gmx::TabulatedNormalDistribution<real, 14> dist;
933 gmx_sd_const_t *sdc = sd->sdc;
934 gmx_sd_sigma_t *sig = sd->sdsig;
936 for (int n = start; n < nrend; n++)
938 int globalAtomIndex = gatindex ? gatindex[n] : n;
939 rng.restart(step, globalAtomIndex);
942 real inverseMass = invmass[n];
943 real invsqrtMass = std::sqrt(inverseMass);
945 int freezeGroup = cFREEZE ? cFREEZE[n] : 0;
946 int accelerationGroup = cACC ? cACC[n] : 0;
947 int temperatureGroup = cTC ? cTC[n] : 0;
949 for (int d = 0; d < DIM; d++)
951 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[freezeGroup][d])
953 if (updateType == SDUpdate::ForcesOnly)
955 real vn = v[n][d] + (inverseMass*f[n][d] + accel[accelerationGroup][d])*dt;
957 // Simple position update.
958 xprime[n][d] = x[n][d] + v[n][d]*dt;
960 else if (updateType == SDUpdate::FrictionAndNoiseOnly)
963 v[n][d] = (vn*sdc[temperatureGroup].em +
964 invsqrtMass*sig[temperatureGroup].V*dist(rng));
965 // The previous phase already updated the
966 // positions with a full v*dt term that must
967 // now be half removed.
968 xprime[n][d] = xprime[n][d] + 0.5*(v[n][d] - vn)*dt;
972 real vn = v[n][d] + (inverseMass*f[n][d] + accel[accelerationGroup][d])*dt;
973 v[n][d] = (vn*sdc[temperatureGroup].em +
974 invsqrtMass*sig[temperatureGroup].V*dist(rng));
975 // Here we include half of the friction+noise
976 // update of v into the position update.
977 xprime[n][d] = x[n][d] + 0.5*(vn + v[n][d])*dt;
982 // When using constraints, the update is split into
983 // two phases, but we only need to zero the update of
984 // virtual, shell or frozen particles in at most one
986 if (updateType != SDUpdate::FrictionAndNoiseOnly)
989 xprime[n][d] = x[n][d];
996 static void do_update_bd(int start, int nrend, real dt,
997 const ivec nFreeze[],
998 const real invmass[], const unsigned short ptype[],
999 const unsigned short cFREEZE[], const unsigned short cTC[],
1000 const rvec x[], rvec xprime[], rvec v[],
1001 const rvec f[], real friction_coefficient,
1002 const real *rf, int64_t step, int seed,
1003 const int* gatindex)
1005 /* note -- these appear to be full step velocities . . . */
1010 // Use 1 bit of internal counters to give us 2*2 64-bits values per stream
1011 // Each 64-bit value is enough for 4 normal distribution table numbers.
1012 gmx::ThreeFry2x64<0> rng(seed, gmx::RandomDomain::UpdateCoordinates);
1013 gmx::TabulatedNormalDistribution<real, 14> dist;
1015 if (friction_coefficient != 0)
1017 invfr = 1.0/friction_coefficient;
1020 for (n = start; (n < nrend); n++)
1022 int ng = gatindex ? gatindex[n] : n;
1024 rng.restart(step, ng);
1035 for (d = 0; (d < DIM); d++)
1037 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
1039 if (friction_coefficient != 0)
1041 vn = invfr*f[n][d] + rf[gt]*dist(rng);
1045 /* NOTE: invmass = 2/(mass*friction_constant*dt) */
1046 vn = 0.5*invmass[n]*f[n][d]*dt
1047 + std::sqrt(0.5*invmass[n])*rf[gt]*dist(rng);
1051 xprime[n][d] = x[n][d]+vn*dt;
1056 xprime[n][d] = x[n][d];
1062 static void calc_ke_part_normal(const rvec v[], const t_grpopts *opts, const t_mdatoms *md,
1063 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel)
1066 t_grp_tcstat *tcstat = ekind->tcstat;
1067 t_grp_acc *grpstat = ekind->grpstat;
1068 int nthread, thread;
1070 /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin. Leap with AveVel is also
1071 an option, but not supported now.
1072 bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
1075 /* group velocities are calculated in update_ekindata and
1076 * accumulated in acumulate_groups.
1077 * Now the partial global and groups ekin.
1079 for (g = 0; (g < opts->ngtc); g++)
1081 copy_mat(tcstat[g].ekinh, tcstat[g].ekinh_old);
1084 clear_mat(tcstat[g].ekinf);
1085 tcstat[g].ekinscalef_nhc = 1.0; /* need to clear this -- logic is complicated! */
1089 clear_mat(tcstat[g].ekinh);
1092 ekind->dekindl_old = ekind->dekindl;
1093 nthread = gmx_omp_nthreads_get(emntUpdate);
1095 #pragma omp parallel for num_threads(nthread) schedule(static)
1096 for (thread = 0; thread < nthread; thread++)
1098 // This OpenMP only loops over arrays and does not call any functions
1099 // or memory allocation. It should not be able to throw, so for now
1100 // we do not need a try/catch wrapper.
1101 int start_t, end_t, n;
1109 start_t = ((thread+0)*md->homenr)/nthread;
1110 end_t = ((thread+1)*md->homenr)/nthread;
1112 ekin_sum = ekind->ekin_work[thread];
1113 dekindl_sum = ekind->dekindl_work[thread];
1115 for (gt = 0; gt < opts->ngtc; gt++)
1117 clear_mat(ekin_sum[gt]);
1123 for (n = start_t; n < end_t; n++)
1133 hm = 0.5*md->massT[n];
1135 for (d = 0; (d < DIM); d++)
1137 v_corrt[d] = v[n][d] - grpstat[ga].u[d];
1139 for (d = 0; (d < DIM); d++)
1141 for (m = 0; (m < DIM); m++)
1143 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1144 ekin_sum[gt][m][d] += hm*v_corrt[m]*v_corrt[d];
1147 if (md->nMassPerturbed && md->bPerturbed[n])
1150 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1156 for (thread = 0; thread < nthread; thread++)
1158 for (g = 0; g < opts->ngtc; g++)
1162 m_add(tcstat[g].ekinf, ekind->ekin_work[thread][g],
1167 m_add(tcstat[g].ekinh, ekind->ekin_work[thread][g],
1172 ekind->dekindl += *ekind->dekindl_work[thread];
1175 inc_nrnb(nrnb, eNR_EKIN, md->homenr);
1178 static void calc_ke_part_visc(const matrix box, const rvec x[], const rvec v[],
1179 const t_grpopts *opts, const t_mdatoms *md,
1180 gmx_ekindata_t *ekind,
1181 t_nrnb *nrnb, gmx_bool bEkinAveVel)
1183 int start = 0, homenr = md->homenr;
1184 int g, d, n, m, gt = 0;
1187 t_grp_tcstat *tcstat = ekind->tcstat;
1188 t_cos_acc *cosacc = &(ekind->cosacc);
1193 for (g = 0; g < opts->ngtc; g++)
1195 copy_mat(ekind->tcstat[g].ekinh, ekind->tcstat[g].ekinh_old);
1196 clear_mat(ekind->tcstat[g].ekinh);
1198 ekind->dekindl_old = ekind->dekindl;
1200 fac = 2*M_PI/box[ZZ][ZZ];
1203 for (n = start; n < start+homenr; n++)
1209 hm = 0.5*md->massT[n];
1211 /* Note that the times of x and v differ by half a step */
1212 /* MRS -- would have to be changed for VV */
1213 cosz = std::cos(fac*x[n][ZZ]);
1214 /* Calculate the amplitude of the new velocity profile */
1215 mvcos += 2*cosz*md->massT[n]*v[n][XX];
1217 copy_rvec(v[n], v_corrt);
1218 /* Subtract the profile for the kinetic energy */
1219 v_corrt[XX] -= cosz*cosacc->vcos;
1220 for (d = 0; (d < DIM); d++)
1222 for (m = 0; (m < DIM); m++)
1224 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1227 tcstat[gt].ekinf[m][d] += hm*v_corrt[m]*v_corrt[d];
1231 tcstat[gt].ekinh[m][d] += hm*v_corrt[m]*v_corrt[d];
1235 if (md->nPerturbed && md->bPerturbed[n])
1237 /* The minus sign here might be confusing.
1238 * The kinetic contribution from dH/dl doesn't come from
1239 * d m(l)/2 v^2 / dl, but rather from d p^2/2m(l) / dl,
1240 * where p are the momenta. The difference is only a minus sign.
1242 dekindl -= 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1245 ekind->dekindl = dekindl;
1246 cosacc->mvcos = mvcos;
1248 inc_nrnb(nrnb, eNR_EKIN, homenr);
1251 void calc_ke_part(const t_state *state, const t_grpopts *opts, const t_mdatoms *md,
1252 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel)
1254 if (ekind->cosacc.cos_accel == 0)
1256 calc_ke_part_normal(state->v.rvec_array(), opts, md, ekind, nrnb, bEkinAveVel);
1260 calc_ke_part_visc(state->box, state->x.rvec_array(), state->v.rvec_array(), opts, md, ekind, nrnb, bEkinAveVel);
1264 extern void init_ekinstate(ekinstate_t *ekinstate, const t_inputrec *ir)
1266 ekinstate->ekin_n = ir->opts.ngtc;
1267 snew(ekinstate->ekinh, ekinstate->ekin_n);
1268 snew(ekinstate->ekinf, ekinstate->ekin_n);
1269 snew(ekinstate->ekinh_old, ekinstate->ekin_n);
1270 ekinstate->ekinscalef_nhc.resize(ekinstate->ekin_n);
1271 ekinstate->ekinscaleh_nhc.resize(ekinstate->ekin_n);
1272 ekinstate->vscale_nhc.resize(ekinstate->ekin_n);
1273 ekinstate->dekindl = 0;
1274 ekinstate->mvcos = 0;
1277 void update_ekinstate(ekinstate_t *ekinstate, const gmx_ekindata_t *ekind)
1281 for (i = 0; i < ekinstate->ekin_n; i++)
1283 copy_mat(ekind->tcstat[i].ekinh, ekinstate->ekinh[i]);
1284 copy_mat(ekind->tcstat[i].ekinf, ekinstate->ekinf[i]);
1285 copy_mat(ekind->tcstat[i].ekinh_old, ekinstate->ekinh_old[i]);
1286 ekinstate->ekinscalef_nhc[i] = ekind->tcstat[i].ekinscalef_nhc;
1287 ekinstate->ekinscaleh_nhc[i] = ekind->tcstat[i].ekinscaleh_nhc;
1288 ekinstate->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
1291 copy_mat(ekind->ekin, ekinstate->ekin_total);
1292 ekinstate->dekindl = ekind->dekindl;
1293 ekinstate->mvcos = ekind->cosacc.mvcos;
1297 void restore_ekinstate_from_state(const t_commrec *cr,
1298 gmx_ekindata_t *ekind, const ekinstate_t *ekinstate)
1304 for (i = 0; i < ekinstate->ekin_n; i++)
1306 copy_mat(ekinstate->ekinh[i], ekind->tcstat[i].ekinh);
1307 copy_mat(ekinstate->ekinf[i], ekind->tcstat[i].ekinf);
1308 copy_mat(ekinstate->ekinh_old[i], ekind->tcstat[i].ekinh_old);
1309 ekind->tcstat[i].ekinscalef_nhc = ekinstate->ekinscalef_nhc[i];
1310 ekind->tcstat[i].ekinscaleh_nhc = ekinstate->ekinscaleh_nhc[i];
1311 ekind->tcstat[i].vscale_nhc = ekinstate->vscale_nhc[i];
1314 copy_mat(ekinstate->ekin_total, ekind->ekin);
1316 ekind->dekindl = ekinstate->dekindl;
1317 ekind->cosacc.mvcos = ekinstate->mvcos;
1318 n = ekinstate->ekin_n;
1323 gmx_bcast(sizeof(n), &n, cr);
1324 for (i = 0; i < n; i++)
1326 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh[0][0]),
1327 ekind->tcstat[i].ekinh[0], cr);
1328 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinf[0][0]),
1329 ekind->tcstat[i].ekinf[0], cr);
1330 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh_old[0][0]),
1331 ekind->tcstat[i].ekinh_old[0], cr);
1333 gmx_bcast(sizeof(ekind->tcstat[i].ekinscalef_nhc),
1334 &(ekind->tcstat[i].ekinscalef_nhc), cr);
1335 gmx_bcast(sizeof(ekind->tcstat[i].ekinscaleh_nhc),
1336 &(ekind->tcstat[i].ekinscaleh_nhc), cr);
1337 gmx_bcast(sizeof(ekind->tcstat[i].vscale_nhc),
1338 &(ekind->tcstat[i].vscale_nhc), cr);
1340 gmx_bcast(DIM*DIM*sizeof(ekind->ekin[0][0]),
1341 ekind->ekin[0], cr);
1343 gmx_bcast(sizeof(ekind->dekindl), &ekind->dekindl, cr);
1344 gmx_bcast(sizeof(ekind->cosacc.mvcos), &ekind->cosacc.mvcos, cr);
1348 void update_tcouple(int64_t step,
1349 const t_inputrec *inputrec,
1351 gmx_ekindata_t *ekind,
1352 const t_extmass *MassQ,
1353 const t_mdatoms *md)
1356 bool doTemperatureCoupling = false;
1358 /* if using vv with trotter decomposition methods, we do this elsewhere in the code */
1359 if (!(EI_VV(inputrec->eI) &&
1360 (inputrecNvtTrotter(inputrec) || inputrecNptTrotter(inputrec) || inputrecNphTrotter(inputrec))))
1362 doTemperatureCoupling = isTemperatureCouplingStep(step, inputrec);
1365 if (doTemperatureCoupling)
1367 real dttc = inputrec->nsttcouple*inputrec->delta_t;
1369 switch (inputrec->etc)
1374 berendsen_tcoupl(inputrec, ekind, dttc, state->therm_integral);
1377 nosehoover_tcoupl(&(inputrec->opts), ekind, dttc,
1378 state->nosehoover_xi.data(), state->nosehoover_vxi.data(), MassQ);
1381 vrescale_tcoupl(inputrec, step, ekind, dttc,
1382 state->therm_integral.data());
1385 /* rescale in place here */
1386 if (EI_VV(inputrec->eI))
1388 rescale_velocities(ekind, md, 0, md->homenr, state->v.rvec_array());
1393 /* Set the T scaling lambda to 1 to have no scaling */
1394 for (int i = 0; (i < inputrec->opts.ngtc); i++)
1396 ekind->tcstat[i].lambda = 1.0;
1401 /*! \brief Computes the atom range for a thread to operate on, ensured SIMD aligned ranges
1403 * \param[in] numThreads The number of threads to divide atoms over
1404 * \param[in] threadIndex The thread to get the range for
1405 * \param[in] numAtoms The total number of atoms (on this rank)
1406 * \param[out] startAtom The start of the atom range
1407 * \param[out] endAtom The end of the atom range, note that this is in general not a multiple of the SIMD width
1409 static void getThreadAtomRange(int numThreads, int threadIndex, int numAtoms,
1410 int *startAtom, int *endAtom)
1412 #if GMX_HAVE_SIMD_UPDATE
1413 constexpr int blockSize = GMX_SIMD_REAL_WIDTH;
1415 constexpr int blockSize = 1;
1417 int numBlocks = (numAtoms + blockSize - 1)/blockSize;
1419 *startAtom = ((numBlocks* threadIndex )/numThreads)*blockSize;
1420 *endAtom = ((numBlocks*(threadIndex + 1))/numThreads)*blockSize;
1421 if (threadIndex == numThreads - 1)
1423 *endAtom = numAtoms;
1427 void update_pcouple_before_coordinates(FILE *fplog,
1429 const t_inputrec *inputrec,
1431 matrix parrinellorahmanMu,
1435 /* Berendsen P-coupling is completely handled after the coordinate update.
1436 * Trotter P-coupling is handled by separate calls to trotter_update().
1438 if (inputrec->epc == epcPARRINELLORAHMAN &&
1439 isPressureCouplingStep(step, inputrec))
1441 real dtpc = inputrec->nstpcouple*inputrec->delta_t;
1443 parrinellorahman_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev,
1444 state->box, state->box_rel, state->boxv,
1445 M, parrinellorahmanMu, bInitStep);
1449 void constrain_velocities(int64_t step,
1450 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1453 gmx::Constraints *constr,
1465 * APPLY CONSTRAINTS:
1472 /* clear out constraints before applying */
1473 clear_mat(vir_part);
1475 /* Constrain the coordinates upd->xp */
1476 constr->apply(do_log, do_ene,
1478 state->x.rvec_array(), state->v.rvec_array(), state->v.rvec_array(),
1480 state->lambda[efptBONDED], dvdlambda,
1481 nullptr, bCalcVir ? &vir_con : nullptr, ConstraintVariable::Velocities);
1485 m_add(vir_part, vir_con, vir_part);
1490 void constrain_coordinates(int64_t step,
1491 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1495 gmx::Constraints *constr,
1508 /* clear out constraints before applying */
1509 clear_mat(vir_part);
1511 /* Constrain the coordinates upd->xp */
1512 constr->apply(do_log, do_ene,
1514 state->x.rvec_array(), upd->xp.rvec_array(), nullptr,
1516 state->lambda[efptBONDED], dvdlambda,
1517 as_rvec_array(state->v.data()), bCalcVir ? &vir_con : nullptr, ConstraintVariable::Positions);
1521 m_add(vir_part, vir_con, vir_part);
1527 update_sd_second_half(int64_t step,
1528 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1529 const t_inputrec *inputrec, /* input record and box stuff */
1530 const t_mdatoms *md,
1532 const t_commrec *cr,
1534 gmx_wallcycle_t wcycle,
1536 gmx::Constraints *constr,
1544 if (inputrec->eI == eiSD1)
1547 int homenr = md->homenr;
1549 /* Cast delta_t from double to real to make the integrators faster.
1550 * The only reason for having delta_t double is to get accurate values
1551 * for t=delta_t*step when step is larger than float precision.
1552 * For integration dt the accuracy of real suffices, since with
1553 * integral += dt*integrand the increment is nearly always (much) smaller
1554 * than the integral (and the integrand has real precision).
1556 real dt = inputrec->delta_t;
1558 wallcycle_start(wcycle, ewcUPDATE);
1560 nth = gmx_omp_nthreads_get(emntUpdate);
1562 #pragma omp parallel for num_threads(nth) schedule(static)
1563 for (th = 0; th < nth; th++)
1567 int start_th, end_th;
1568 getThreadAtomRange(nth, th, homenr, &start_th, &end_th);
1570 doSDUpdateGeneral<SDUpdate::FrictionAndNoiseOnly>
1572 start_th, end_th, dt,
1573 inputrec->opts.acc, inputrec->opts.nFreeze,
1574 md->invmass, md->ptype,
1575 md->cFREEZE, nullptr, md->cTC,
1576 state->x.rvec_array(), upd->xp.rvec_array(),
1577 state->v.rvec_array(), nullptr,
1578 step, inputrec->ld_seed,
1579 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1581 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1583 inc_nrnb(nrnb, eNR_UPDATE, homenr);
1584 wallcycle_stop(wcycle, ewcUPDATE);
1586 /* Constrain the coordinates upd->xp for half a time step */
1587 constr->apply(do_log, do_ene,
1589 state->x.rvec_array(), upd->xp.rvec_array(), nullptr,
1591 state->lambda[efptBONDED], dvdlambda,
1592 as_rvec_array(state->v.data()), nullptr, ConstraintVariable::Positions);
1596 void finish_update(const t_inputrec *inputrec, /* input record and box stuff */
1597 const t_mdatoms *md,
1599 const t_graph *graph,
1601 gmx_wallcycle_t wcycle,
1603 const gmx::Constraints *constr)
1605 int homenr = md->homenr;
1607 /* We must always unshift after updating coordinates; if we did not shake
1608 x was shifted in do_force */
1610 /* NOTE Currently we always integrate to a temporary buffer and
1611 * then copy the results back. */
1613 wallcycle_start_nocount(wcycle, ewcUPDATE);
1615 if (md->cFREEZE != nullptr && constr != nullptr)
1617 /* If we have atoms that are frozen along some, but not all
1618 * dimensions, then any constraints will have moved them also along
1619 * the frozen dimensions. To freeze such degrees of freedom
1620 * we copy them back here to later copy them forward. It would
1621 * be more elegant and slightly more efficient to copies zero
1622 * times instead of twice, but the graph case below prevents this.
1624 const ivec *nFreeze = inputrec->opts.nFreeze;
1625 bool partialFreezeAndConstraints = false;
1626 for (int g = 0; g < inputrec->opts.ngfrz; g++)
1628 int numFreezeDim = nFreeze[g][XX] + nFreeze[g][YY] + nFreeze[g][ZZ];
1629 if (numFreezeDim > 0 && numFreezeDim < 3)
1631 partialFreezeAndConstraints = true;
1634 if (partialFreezeAndConstraints)
1636 auto xp = makeArrayRef(upd->xp).subArray(0, homenr);
1637 auto x = makeConstArrayRef(state->x).subArray(0, homenr);
1638 for (int i = 0; i < homenr; i++)
1640 int g = md->cFREEZE[i];
1642 for (int d = 0; d < DIM; d++)
1653 if (graph && (graph->nnodes > 0))
1655 unshift_x(graph, state->box, state->x.rvec_array(), upd->xp.rvec_array());
1656 if (TRICLINIC(state->box))
1658 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
1662 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
1667 auto xp = makeConstArrayRef(upd->xp).subArray(0, homenr);
1668 auto x = makeArrayRef(state->x).subArray(0, homenr);
1669 #ifndef __clang_analyzer__
1670 int gmx_unused nth = gmx_omp_nthreads_get(emntUpdate);
1672 #pragma omp parallel for num_threads(nth) schedule(static)
1673 for (int i = 0; i < homenr; i++)
1675 // Trivial statement, does not throw
1679 wallcycle_stop(wcycle, ewcUPDATE);
1681 /* ############# END the update of velocities and positions ######### */
1684 void update_pcouple_after_coordinates(FILE *fplog,
1686 const t_inputrec *inputrec,
1687 const t_mdatoms *md,
1688 const matrix pressure,
1689 const matrix forceVirial,
1690 const matrix constraintVirial,
1691 const matrix parrinellorahmanMu,
1697 int homenr = md->homenr;
1699 /* Cast to real for faster code, no loss in precision (see comment above) */
1700 real dt = inputrec->delta_t;
1703 /* now update boxes */
1704 switch (inputrec->epc)
1708 case (epcBERENDSEN):
1709 if (isPressureCouplingStep(step, inputrec))
1711 real dtpc = inputrec->nstpcouple*dt;
1713 berendsen_pcoupl(fplog, step, inputrec, dtpc,
1714 pressure, state->box,
1715 forceVirial, constraintVirial,
1716 mu, &state->baros_integral);
1717 berendsen_pscale(inputrec, mu, state->box, state->box_rel,
1718 start, homenr, state->x.rvec_array(),
1722 case (epcPARRINELLORAHMAN):
1723 if (isPressureCouplingStep(step, inputrec))
1725 /* The box velocities were updated in do_pr_pcoupl,
1726 * but we dont change the box vectors until we get here
1727 * since we need to be able to shift/unshift above.
1729 real dtpc = inputrec->nstpcouple*dt;
1730 for (int i = 0; i < DIM; i++)
1732 for (int m = 0; m <= i; m++)
1734 state->box[i][m] += dtpc*state->boxv[i][m];
1737 preserve_box_shape(inputrec, state->box_rel, state->box);
1739 /* Scale the coordinates */
1740 auto x = state->x.rvec_array();
1741 for (int n = start; n < start + homenr; n++)
1743 tmvmul_ur0(parrinellorahmanMu, x[n], x[n]);
1748 switch (inputrec->epct)
1750 case (epctISOTROPIC):
1751 /* DIM * eta = ln V. so DIM*eta_new = DIM*eta_old + DIM*dt*veta =>
1752 ln V_new = ln V_old + 3*dt*veta => V_new = V_old*exp(3*dt*veta) =>
1753 Side length scales as exp(veta*dt) */
1755 msmul(state->box, std::exp(state->veta*dt), state->box);
1757 /* Relate veta to boxv. veta = d(eta)/dT = (1/DIM)*1/V dV/dT.
1758 o If we assume isotropic scaling, and box length scaling
1759 factor L, then V = L^DIM (det(M)). So dV/dt = DIM
1760 L^(DIM-1) dL/dt det(M), and veta = (1/L) dL/dt. The
1761 determinant of B is L^DIM det(M), and the determinant
1762 of dB/dt is (dL/dT)^DIM det (M). veta will be
1763 (det(dB/dT)/det(B))^(1/3). Then since M =
1764 B_new*(vol_new)^(1/3), dB/dT_new = (veta_new)*B(new). */
1766 msmul(state->box, state->veta, state->boxv);
1778 auto localX = makeArrayRef(state->x).subArray(start, homenr);
1779 upd->deform->apply(localX, state->box, step);
1783 void update_coords(int64_t step,
1784 const t_inputrec *inputrec, /* input record and box stuff */
1785 const t_mdatoms *md,
1787 gmx::ArrayRefWithPadding<gmx::RVec> f,
1788 const t_fcdata *fcd,
1789 const gmx_ekindata_t *ekind,
1793 const t_commrec *cr, /* these shouldn't be here -- need to think about it */
1794 const gmx::Constraints *constr)
1796 gmx_bool bDoConstr = (nullptr != constr);
1798 /* Running the velocity half does nothing except for velocity verlet */
1799 if ((UpdatePart == etrtVELOCITY1 || UpdatePart == etrtVELOCITY2) &&
1800 !EI_VV(inputrec->eI))
1802 gmx_incons("update_coords called for velocity without VV integrator");
1805 int homenr = md->homenr;
1807 /* Cast to real for faster code, no loss in precision (see comment above) */
1808 real dt = inputrec->delta_t;
1810 /* We need to update the NMR restraint history when time averaging is used */
1811 if (state->flags & (1<<estDISRE_RM3TAV))
1813 update_disres_history(fcd, &state->hist);
1815 if (state->flags & (1<<estORIRE_DTAV))
1817 update_orires_history(fcd, &state->hist);
1820 /* ############# START The update of velocities and positions ######### */
1821 int nth = gmx_omp_nthreads_get(emntUpdate);
1823 #pragma omp parallel for num_threads(nth) schedule(static)
1824 for (int th = 0; th < nth; th++)
1828 int start_th, end_th;
1829 getThreadAtomRange(nth, th, homenr, &start_th, &end_th);
1831 const rvec *x_rvec = state->x.rvec_array();
1832 rvec *xp_rvec = upd->xp.rvec_array();
1833 rvec *v_rvec = state->v.rvec_array();
1834 const rvec *f_rvec = as_rvec_array(f.unpaddedArrayRef().data());
1836 switch (inputrec->eI)
1839 do_update_md(start_th, end_th, step, dt,
1840 inputrec, md, ekind, state->box,
1841 x_rvec, xp_rvec, v_rvec, f_rvec,
1842 state->nosehoover_vxi.data(), M);
1847 // With constraints, the SD update is done in 2 parts
1848 doSDUpdateGeneral<SDUpdate::ForcesOnly>
1850 start_th, end_th, dt,
1851 inputrec->opts.acc, inputrec->opts.nFreeze,
1852 md->invmass, md->ptype,
1853 md->cFREEZE, md->cACC, nullptr,
1854 x_rvec, xp_rvec, v_rvec, f_rvec,
1855 step, inputrec->ld_seed, nullptr);
1859 doSDUpdateGeneral<SDUpdate::Combined>
1861 start_th, end_th, dt,
1862 inputrec->opts.acc, inputrec->opts.nFreeze,
1863 md->invmass, md->ptype,
1864 md->cFREEZE, md->cACC, md->cTC,
1865 x_rvec, xp_rvec, v_rvec, f_rvec,
1866 step, inputrec->ld_seed,
1867 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1871 do_update_bd(start_th, end_th, dt,
1872 inputrec->opts.nFreeze, md->invmass, md->ptype,
1873 md->cFREEZE, md->cTC,
1874 x_rvec, xp_rvec, v_rvec, f_rvec,
1877 step, inputrec->ld_seed, DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1882 gmx_bool bExtended = (inputrec->etc == etcNOSEHOOVER ||
1883 inputrec->epc == epcPARRINELLORAHMAN ||
1884 inputrec->epc == epcMTTK);
1886 /* assuming barostat coupled to group 0 */
1887 real alpha = 1.0 + DIM/static_cast<real>(inputrec->opts.nrdf[0]);
1892 do_update_vv_vel(start_th, end_th, dt,
1893 inputrec->opts.acc, inputrec->opts.nFreeze,
1894 md->invmass, md->ptype,
1895 md->cFREEZE, md->cACC,
1897 bExtended, state->veta, alpha);
1900 do_update_vv_pos(start_th, end_th, dt,
1901 inputrec->opts.nFreeze,
1902 md->ptype, md->cFREEZE,
1903 x_rvec, xp_rvec, v_rvec,
1904 bExtended, state->veta);
1910 gmx_fatal(FARGS, "Don't know how to update coordinates");
1913 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1918 extern gmx_bool update_randomize_velocities(const t_inputrec *ir, int64_t step, const t_commrec *cr,
1919 const t_mdatoms *md,
1920 gmx::ArrayRef<gmx::RVec> v,
1921 const gmx_update_t *upd,
1922 const gmx::Constraints *constr)
1925 real rate = (ir->delta_t)/ir->opts.tau_t[0];
1927 if (ir->etc == etcANDERSEN && constr != nullptr)
1929 /* Currently, Andersen thermostat does not support constrained
1930 systems. Functionality exists in the andersen_tcoupl
1931 function in GROMACS 4.5.7 to allow this combination. That
1932 code could be ported to the current random-number
1933 generation approach, but has not yet been done because of
1934 lack of time and resources. */
1935 gmx_fatal(FARGS, "Normal Andersen is currently not supported with constraints, use massive Andersen instead");
1938 /* proceed with andersen if 1) it's fixed probability per
1939 particle andersen or 2) it's massive andersen and it's tau_t/dt */
1940 if ((ir->etc == etcANDERSEN) || do_per_step(step, roundToInt(1.0/rate)))
1942 andersen_tcoupl(ir, step, cr, md, v, rate,
1943 upd->sd->randomize_group, upd->sd->boltzfac);