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 */
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 rvec accel[], 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,
723 const unsigned short ptype[], const unsigned short cFREEZE[],
724 const rvec x[], rvec xprime[], 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);
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.resize(gmx::paddedRVecVectorSize(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 rvec accel[], 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 if (updateType != SDUpdate::FrictionAndNoiseOnly)
913 GMX_ASSERT(f != nullptr, "SD update with forces requires forces");
914 GMX_ASSERT(cACC != nullptr, "SD update with forces requires acceleration groups");
916 if (updateType != SDUpdate::ForcesOnly)
918 GMX_ASSERT(cTC != nullptr, "SD update with noise requires temperature groups");
921 // Even 0 bits internal counter gives 2x64 ints (more than enough for three table lookups)
922 gmx::ThreeFry2x64<0> rng(seed, gmx::RandomDomain::UpdateCoordinates);
923 gmx::TabulatedNormalDistribution<real, 14> dist;
925 gmx_sd_const_t *sdc = sd->sdc;
926 gmx_sd_sigma_t *sig = sd->sdsig;
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*sdc[temperatureGroup].em +
956 invsqrtMass*sig[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*sdc[temperatureGroup].em +
966 invsqrtMass*sig[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,
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(rvec v[], t_grpopts *opts, 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(matrix box, rvec x[], rvec v[],
1171 t_grpopts *opts, 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(t_state *state, t_grpopts *opts, 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(as_rvec_array(state->v.data()), opts, md, ekind, nrnb, bEkinAveVel);
1252 calc_ke_part_visc(state->box, as_rvec_array(state->x.data()), as_rvec_array(state->v.data()), 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, 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 t_inputrec *inputrec,
1343 gmx_ekindata_t *ekind,
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, as_rvec_array(state->v.data()));
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_wallcycle_t wcycle,
1446 gmx::Constraints *constr,
1458 * APPLY CONSTRAINTS:
1465 /* clear out constraints before applying */
1466 clear_mat(vir_part);
1468 /* Constrain the coordinates upd->xp */
1469 wallcycle_start(wcycle, ewcCONSTR);
1471 constr->apply(do_log, do_ene,
1473 as_rvec_array(state->x.data()), as_rvec_array(state->v.data()), as_rvec_array(state->v.data()),
1475 state->lambda[efptBONDED], dvdlambda,
1476 nullptr, bCalcVir ? &vir_con : nullptr, ConstraintVariable::Velocities);
1478 wallcycle_stop(wcycle, ewcCONSTR);
1482 m_add(vir_part, vir_con, vir_part);
1487 void constrain_coordinates(int64_t step,
1488 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1491 gmx_wallcycle_t wcycle,
1493 gmx::Constraints *constr,
1506 /* clear out constraints before applying */
1507 clear_mat(vir_part);
1509 /* Constrain the coordinates upd->xp */
1510 wallcycle_start(wcycle, ewcCONSTR);
1512 constr->apply(do_log, do_ene,
1514 as_rvec_array(state->x.data()), as_rvec_array(upd->xp.data()), nullptr,
1516 state->lambda[efptBONDED], dvdlambda,
1517 as_rvec_array(state->v.data()), bCalcVir ? &vir_con : nullptr, ConstraintVariable::Positions);
1519 wallcycle_stop(wcycle, ewcCONSTR);
1523 m_add(vir_part, vir_con, vir_part);
1529 update_sd_second_half(int64_t step,
1530 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1531 const t_inputrec *inputrec, /* input record and box stuff */
1534 const t_commrec *cr,
1536 gmx_wallcycle_t wcycle,
1538 gmx::Constraints *constr,
1546 if (inputrec->eI == eiSD1)
1549 int homenr = md->homenr;
1551 /* Cast delta_t from double to real to make the integrators faster.
1552 * The only reason for having delta_t double is to get accurate values
1553 * for t=delta_t*step when step is larger than float precision.
1554 * For integration dt the accuracy of real suffices, since with
1555 * integral += dt*integrand the increment is nearly always (much) smaller
1556 * than the integral (and the integrand has real precision).
1558 real dt = inputrec->delta_t;
1560 wallcycle_start(wcycle, ewcUPDATE);
1562 nth = gmx_omp_nthreads_get(emntUpdate);
1564 #pragma omp parallel for num_threads(nth) schedule(static)
1565 for (th = 0; th < nth; th++)
1569 int start_th, end_th;
1570 getThreadAtomRange(nth, th, homenr, &start_th, &end_th);
1572 doSDUpdateGeneral<SDUpdate::FrictionAndNoiseOnly>
1574 start_th, end_th, dt,
1575 inputrec->opts.acc, inputrec->opts.nFreeze,
1576 md->invmass, md->ptype,
1577 md->cFREEZE, nullptr, md->cTC,
1578 as_rvec_array(state->x.data()), as_rvec_array(upd->xp.data()),
1579 as_rvec_array(state->v.data()), nullptr,
1580 step, inputrec->ld_seed,
1581 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1583 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1585 inc_nrnb(nrnb, eNR_UPDATE, homenr);
1586 wallcycle_stop(wcycle, ewcUPDATE);
1589 /* Constrain the coordinates upd->xp for half a time step */
1590 wallcycle_start(wcycle, ewcCONSTR);
1591 constr->apply(do_log, do_ene,
1593 as_rvec_array(state->x.data()), as_rvec_array(upd->xp.data()), nullptr,
1595 state->lambda[efptBONDED], dvdlambda,
1596 as_rvec_array(state->v.data()), nullptr, ConstraintVariable::Positions);
1597 wallcycle_stop(wcycle, ewcCONSTR);
1602 void finish_update(const t_inputrec *inputrec, /* input record and box stuff */
1607 gmx_wallcycle_t wcycle,
1609 gmx::Constraints *constr)
1611 int homenr = md->homenr;
1613 /* We must always unshift after updating coordinates; if we did not shake
1614 x was shifted in do_force */
1616 /* NOTE Currently we always integrate to a temporary buffer and
1617 * then copy the results back. */
1619 wallcycle_start_nocount(wcycle, ewcUPDATE);
1621 if (md->cFREEZE != nullptr && constr != nullptr)
1623 /* If we have atoms that are frozen along some, but not all
1624 * dimensions, then any constraints will have moved them also along
1625 * the frozen dimensions. To freeze such degrees of freedom
1626 * we copy them back here to later copy them forward. It would
1627 * be more elegant and slightly more efficient to copies zero
1628 * times instead of twice, but the graph case below prevents this.
1630 const ivec *nFreeze = inputrec->opts.nFreeze;
1631 bool partialFreezeAndConstraints = false;
1632 for (int g = 0; g < inputrec->opts.ngfrz; g++)
1634 int numFreezeDim = nFreeze[g][XX] + nFreeze[g][YY] + nFreeze[g][ZZ];
1635 if (numFreezeDim > 0 && numFreezeDim < 3)
1637 partialFreezeAndConstraints = true;
1640 if (partialFreezeAndConstraints)
1642 for (int i = 0; i < homenr; i++)
1644 int g = md->cFREEZE[i];
1646 for (int d = 0; d < DIM; d++)
1650 upd->xp[i][d] = state->x[i][d];
1657 if (graph && (graph->nnodes > 0))
1659 unshift_x(graph, state->box, as_rvec_array(state->x.data()), as_rvec_array(upd->xp.data()));
1660 if (TRICLINIC(state->box))
1662 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
1666 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
1671 /* The copy is performance sensitive, so use a bare pointer */
1672 rvec *xp = as_rvec_array(upd->xp.data());
1673 #ifndef __clang_analyzer__
1674 int gmx_unused nth = gmx_omp_nthreads_get(emntUpdate);
1676 #pragma omp parallel for num_threads(nth) schedule(static)
1677 for (int i = 0; i < homenr; i++)
1679 // Trivial statement, does not throw
1680 copy_rvec(xp[i], state->x[i]);
1683 wallcycle_stop(wcycle, ewcUPDATE);
1685 /* ############# END the update of velocities and positions ######### */
1688 void update_pcouple_after_coordinates(FILE *fplog,
1690 const t_inputrec *inputrec,
1691 const t_mdatoms *md,
1692 const matrix pressure,
1693 const matrix forceVirial,
1694 const matrix constraintVirial,
1695 const matrix parrinellorahmanMu,
1701 int homenr = md->homenr;
1703 /* Cast to real for faster code, no loss in precision (see comment above) */
1704 real dt = inputrec->delta_t;
1707 /* now update boxes */
1708 switch (inputrec->epc)
1712 case (epcBERENDSEN):
1713 if (isPressureCouplingStep(step, inputrec))
1715 real dtpc = inputrec->nstpcouple*dt;
1717 berendsen_pcoupl(fplog, step, inputrec, dtpc,
1718 pressure, state->box,
1719 forceVirial, constraintVirial,
1720 mu, &state->baros_integral);
1721 berendsen_pscale(inputrec, mu, state->box, state->box_rel,
1722 start, homenr, as_rvec_array(state->x.data()),
1726 case (epcPARRINELLORAHMAN):
1727 if (isPressureCouplingStep(step, inputrec))
1729 /* The box velocities were updated in do_pr_pcoupl,
1730 * but we dont change the box vectors until we get here
1731 * since we need to be able to shift/unshift above.
1733 real dtpc = inputrec->nstpcouple*dt;
1734 for (int i = 0; i < DIM; i++)
1736 for (int m = 0; m <= i; m++)
1738 state->box[i][m] += dtpc*state->boxv[i][m];
1741 preserve_box_shape(inputrec, state->box_rel, state->box);
1743 /* Scale the coordinates */
1744 for (int n = start; n < start + homenr; n++)
1746 tmvmul_ur0(parrinellorahmanMu, state->x[n], state->x[n]);
1751 switch (inputrec->epct)
1753 case (epctISOTROPIC):
1754 /* DIM * eta = ln V. so DIM*eta_new = DIM*eta_old + DIM*dt*veta =>
1755 ln V_new = ln V_old + 3*dt*veta => V_new = V_old*exp(3*dt*veta) =>
1756 Side length scales as exp(veta*dt) */
1758 msmul(state->box, std::exp(state->veta*dt), state->box);
1760 /* Relate veta to boxv. veta = d(eta)/dT = (1/DIM)*1/V dV/dT.
1761 o If we assume isotropic scaling, and box length scaling
1762 factor L, then V = L^DIM (det(M)). So dV/dt = DIM
1763 L^(DIM-1) dL/dt det(M), and veta = (1/L) dL/dt. The
1764 determinant of B is L^DIM det(M), and the determinant
1765 of dB/dt is (dL/dT)^DIM det (M). veta will be
1766 (det(dB/dT)/det(B))^(1/3). Then since M =
1767 B_new*(vol_new)^(1/3), dB/dT_new = (veta_new)*B(new). */
1769 msmul(state->box, state->veta, state->boxv);
1781 auto localX = gmx::ArrayRef<gmx::RVec>(state->x).subArray(start, homenr);
1782 upd->deform->apply(localX, state->box, step);
1786 void update_coords(int64_t step,
1787 t_inputrec *inputrec, /* input record and box stuff */
1790 gmx::PaddedArrayRef<gmx::RVec> f, /* forces on home particles */
1792 gmx_ekindata_t *ekind,
1796 const t_commrec *cr, /* these shouldn't be here -- need to think about it */
1797 gmx::Constraints *constr)
1799 gmx_bool bDoConstr = (nullptr != constr);
1801 /* Running the velocity half does nothing except for velocity verlet */
1802 if ((UpdatePart == etrtVELOCITY1 || UpdatePart == etrtVELOCITY2) &&
1803 !EI_VV(inputrec->eI))
1805 gmx_incons("update_coords called for velocity without VV integrator");
1808 int homenr = md->homenr;
1810 /* Cast to real for faster code, no loss in precision (see comment above) */
1811 real dt = inputrec->delta_t;
1813 /* We need to update the NMR restraint history when time averaging is used */
1814 if (state->flags & (1<<estDISRE_RM3TAV))
1816 update_disres_history(fcd, &state->hist);
1818 if (state->flags & (1<<estORIRE_DTAV))
1820 update_orires_history(fcd, &state->hist);
1823 /* ############# START The update of velocities and positions ######### */
1824 int nth = gmx_omp_nthreads_get(emntUpdate);
1826 #pragma omp parallel for num_threads(nth) schedule(static)
1827 for (int th = 0; th < nth; th++)
1831 int start_th, end_th;
1832 getThreadAtomRange(nth, th, homenr, &start_th, &end_th);
1834 /* Strictly speaking, we would only need this check with SIMD
1835 * and for the actual SIMD width. But since the code currently
1836 * always adds padding for GMX_REAL_MAX_SIMD_WIDTH, we check that.
1838 gmx::index gmx_used_in_debug homenrSimdPadded;
1839 homenrSimdPadded = ((homenr + GMX_REAL_MAX_SIMD_WIDTH - 1)/GMX_REAL_MAX_SIMD_WIDTH)*GMX_REAL_MAX_SIMD_WIDTH;
1840 GMX_ASSERT(gmx::index(state->x.size()) >= homenrSimdPadded, "state->x needs to be padded for SIMD access");
1841 GMX_ASSERT(gmx::index(upd->xp.size()) >= homenrSimdPadded, "upd->xp needs to be padded for SIMD access");
1842 GMX_ASSERT(gmx::index(state->v.size()) >= homenrSimdPadded, "state->v needs to be padded for SIMD access");
1843 GMX_ASSERT(f.size() >= homenrSimdPadded, "f needs to be padded for SIMD access");
1845 const rvec *x_rvec = as_rvec_array(state->x.data());
1846 rvec *xp_rvec = as_rvec_array(upd->xp.data());
1847 rvec *v_rvec = as_rvec_array(state->v.data());
1848 const rvec *f_rvec = as_rvec_array(f.data());
1850 switch (inputrec->eI)
1853 do_update_md(start_th, end_th, step, dt,
1854 inputrec, md, ekind, state->box,
1855 x_rvec, xp_rvec, v_rvec, f_rvec,
1856 state->nosehoover_vxi.data(), M);
1861 // With constraints, the SD update is done in 2 parts
1862 doSDUpdateGeneral<SDUpdate::ForcesOnly>
1864 start_th, end_th, dt,
1865 inputrec->opts.acc, inputrec->opts.nFreeze,
1866 md->invmass, md->ptype,
1867 md->cFREEZE, md->cACC, nullptr,
1868 x_rvec, xp_rvec, v_rvec, f_rvec,
1869 step, inputrec->ld_seed, nullptr);
1873 doSDUpdateGeneral<SDUpdate::Combined>
1875 start_th, end_th, dt,
1876 inputrec->opts.acc, inputrec->opts.nFreeze,
1877 md->invmass, md->ptype,
1878 md->cFREEZE, md->cACC, md->cTC,
1879 x_rvec, xp_rvec, v_rvec, f_rvec,
1880 step, inputrec->ld_seed,
1881 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1885 do_update_bd(start_th, end_th, dt,
1886 inputrec->opts.nFreeze, md->invmass, md->ptype,
1887 md->cFREEZE, md->cTC,
1888 x_rvec, xp_rvec, v_rvec, f_rvec,
1891 step, inputrec->ld_seed, DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1896 gmx_bool bExtended = (inputrec->etc == etcNOSEHOOVER ||
1897 inputrec->epc == epcPARRINELLORAHMAN ||
1898 inputrec->epc == epcMTTK);
1900 /* assuming barostat coupled to group 0 */
1901 real alpha = 1.0 + DIM/static_cast<real>(inputrec->opts.nrdf[0]);
1906 do_update_vv_vel(start_th, end_th, dt,
1907 inputrec->opts.acc, inputrec->opts.nFreeze,
1908 md->invmass, md->ptype,
1909 md->cFREEZE, md->cACC,
1911 bExtended, state->veta, alpha);
1914 do_update_vv_pos(start_th, end_th, dt,
1915 inputrec->opts.nFreeze,
1916 md->ptype, md->cFREEZE,
1917 x_rvec, xp_rvec, v_rvec,
1918 bExtended, state->veta);
1924 gmx_fatal(FARGS, "Don't know how to update coordinates");
1927 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1932 extern gmx_bool update_randomize_velocities(t_inputrec *ir, int64_t step, const t_commrec *cr,
1933 t_mdatoms *md, t_state *state, gmx_update_t *upd, gmx::Constraints *constr)
1936 real rate = (ir->delta_t)/ir->opts.tau_t[0];
1938 if (ir->etc == etcANDERSEN && constr != nullptr)
1940 /* Currently, Andersen thermostat does not support constrained
1941 systems. Functionality exists in the andersen_tcoupl
1942 function in GROMACS 4.5.7 to allow this combination. That
1943 code could be ported to the current random-number
1944 generation approach, but has not yet been done because of
1945 lack of time and resources. */
1946 gmx_fatal(FARGS, "Normal Andersen is currently not supported with constraints, use massive Andersen instead");
1949 /* proceed with andersen if 1) it's fixed probability per
1950 particle andersen or 2) it's massive andersen and it's tau_t/dt */
1951 if ((ir->etc == etcANDERSEN) || do_per_step(step, static_cast<int>(1.0/rate + 0.5)))
1953 andersen_tcoupl(ir, step, cr, md, state, rate,
1954 upd->sd->randomize_group, upd->sd->boltzfac);