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/compat/make_unique.h"
47 #include "gromacs/domdec/domdec_struct.h"
48 #include "gromacs/fileio/confio.h"
49 #include "gromacs/gmxlib/network.h"
50 #include "gromacs/gmxlib/nrnb.h"
51 #include "gromacs/listed-forces/disre.h"
52 #include "gromacs/listed-forces/orires.h"
53 #include "gromacs/math/functions.h"
54 #include "gromacs/math/invertmatrix.h"
55 #include "gromacs/math/paddedvector.h"
56 #include "gromacs/math/units.h"
57 #include "gromacs/math/vec.h"
58 #include "gromacs/math/vecdump.h"
59 #include "gromacs/mdlib/boxdeformation.h"
60 #include "gromacs/mdlib/constr.h"
61 #include "gromacs/mdlib/gmx_omp_nthreads.h"
62 #include "gromacs/mdlib/mdrun.h"
63 #include "gromacs/mdlib/sim_util.h"
64 #include "gromacs/mdlib/tgroup.h"
65 #include "gromacs/mdtypes/commrec.h"
66 #include "gromacs/mdtypes/group.h"
67 #include "gromacs/mdtypes/inputrec.h"
68 #include "gromacs/mdtypes/md_enums.h"
69 #include "gromacs/mdtypes/state.h"
70 #include "gromacs/pbcutil/boxutilities.h"
71 #include "gromacs/pbcutil/mshift.h"
72 #include "gromacs/pbcutil/pbc.h"
73 #include "gromacs/pulling/pull.h"
74 #include "gromacs/random/tabulatednormaldistribution.h"
75 #include "gromacs/random/threefry.h"
76 #include "gromacs/simd/simd.h"
77 #include "gromacs/timing/wallcycle.h"
78 #include "gromacs/topology/atoms.h"
79 #include "gromacs/utility/exceptions.h"
80 #include "gromacs/utility/fatalerror.h"
81 #include "gromacs/utility/futil.h"
82 #include "gromacs/utility/gmxassert.h"
83 #include "gromacs/utility/gmxomp.h"
84 #include "gromacs/utility/smalloc.h"
86 using namespace gmx; // TODO: Remove when this file is moved into gmx namespace
99 std::vector<real> bd_rf;
101 std::vector<gmx_sd_const_t> sdc;
102 std::vector<gmx_sd_sigma_t> sdsig;
103 /* andersen temperature control stuff */
104 std::vector<bool> randomize_group;
105 std::vector<real> boltzfac;
107 gmx_stochd_t(const t_inputrec *ir);
112 std::unique_ptr<gmx_stochd_t> sd;
113 /* xprime for constraint algorithms */
114 PaddedVector<gmx::RVec> xp;
116 /* Variables for the deform algorithm */
117 int64_t deformref_step;
118 matrix deformref_box;
120 //! Box deformation handler (or nullptr if inactive).
121 gmx::BoxDeformation *deform;
124 static bool isTemperatureCouplingStep(int64_t step, const t_inputrec *ir)
126 /* We should only couple after a step where energies were determined (for leapfrog versions)
127 or the step energies are determined, for velocity verlet versions */
137 return ir->etc != etcNO &&
138 (ir->nsttcouple == 1 ||
139 do_per_step(step + ir->nsttcouple - offset, ir->nsttcouple));
142 static bool isPressureCouplingStep(int64_t step, const t_inputrec *ir)
144 GMX_ASSERT(ir->epc != epcMTTK, "MTTK pressure coupling is not handled here");
147 if (ir->epc == epcBERENDSEN)
155 /* We should only couple after a step where pressures were determined */
156 return ir->epc != etcNO &&
157 (ir->nstpcouple == 1 ||
158 do_per_step(step + ir->nstpcouple - offset, ir->nstpcouple));
161 /*! \brief Sets the velocities of virtual sites to zero */
162 static void clearVsiteVelocities(int start,
164 const unsigned short *particleType,
165 rvec * gmx_restrict v)
167 for (int a = start; a < nrend; a++)
169 if (particleType[a] == eptVSite)
176 /*! \brief Sets the number of different temperature coupling values */
177 enum class NumTempScaleValues
179 single, //!< Single T-scaling value (either one group or all values =1)
180 multiple //!< Multiple T-scaling values, need to use T-group indices
183 /*! \brief Sets if to apply no or diagonal Parrinello-Rahman pressure scaling
185 * Note that this enum is only used in updateMDLeapfrogSimple(), which does
186 * not handle fully anistropic Parrinello-Rahman scaling, so we only have
187 * options \p no and \p diagonal here and no anistropic option.
189 enum class ApplyParrinelloRahmanVScaling
191 no, //!< Do not apply velocity scaling (not a PR-coupling run or step)
192 diagonal //!< Apply velocity scaling using a diagonal matrix
195 /*! \brief Integrate using leap-frog with T-scaling and optionally diagonal Parrinello-Rahman p-coupling
197 * \tparam numTempScaleValues The number of different T-couple values
198 * \tparam applyPRVScaling Apply Parrinello-Rahman velocity scaling
199 * \param[in] start Index of first atom to update
200 * \param[in] nrend Last atom to update: \p nrend - 1
201 * \param[in] dt The time step
202 * \param[in] dtPressureCouple Time step for pressure coupling
203 * \param[in] invMassPerDim 1/mass per atom and dimension
204 * \param[in] tcstat Temperature coupling information
205 * \param[in] cTC T-coupling group index per atom
206 * \param[in] pRVScaleMatrixDiagonal Parrinello-Rahman v-scale matrix diagonal
207 * \param[in] x Input coordinates
208 * \param[out] xprime Updated coordinates
209 * \param[inout] v Velocities
210 * \param[in] f Forces
212 * We expect this template to get good SIMD acceleration by most compilers,
213 * unlike the more complex general template.
214 * Note that we might get even better SIMD acceleration when we introduce
215 * aligned (and padded) memory, possibly with some hints for the compilers.
217 template<NumTempScaleValues numTempScaleValues,
218 ApplyParrinelloRahmanVScaling applyPRVScaling>
220 updateMDLeapfrogSimple(int start,
223 real dtPressureCouple,
224 const rvec * gmx_restrict invMassPerDim,
225 const t_grp_tcstat * tcstat,
226 const unsigned short * cTC,
227 const rvec pRVScaleMatrixDiagonal,
228 const rvec * gmx_restrict x,
229 rvec * gmx_restrict xprime,
230 rvec * gmx_restrict v,
231 const rvec * gmx_restrict f)
235 if (numTempScaleValues == NumTempScaleValues::single)
237 lambdaGroup = tcstat[0].lambda;
240 for (int a = start; a < nrend; a++)
242 if (numTempScaleValues == NumTempScaleValues::multiple)
244 lambdaGroup = tcstat[cTC[a]].lambda;
247 for (int d = 0; d < DIM; d++)
249 /* Note that using rvec invMassPerDim results in more efficient
250 * SIMD code, but this increases the cache pressure.
251 * For large systems with PME on the CPU this slows down the
252 * (then already slow) update by 20%. If all data remains in cache,
253 * using rvec is much faster.
255 real vNew = lambdaGroup*v[a][d] + f[a][d]*invMassPerDim[a][d]*dt;
257 if (applyPRVScaling == ApplyParrinelloRahmanVScaling::diagonal)
259 vNew -= dtPressureCouple*pRVScaleMatrixDiagonal[d]*v[a][d];
262 xprime[a][d] = x[a][d] + vNew*dt;
267 #if GMX_SIMD && GMX_SIMD_HAVE_REAL
268 #define GMX_HAVE_SIMD_UPDATE 1
270 #define GMX_HAVE_SIMD_UPDATE 0
273 #if GMX_HAVE_SIMD_UPDATE
275 /*! \brief Load (aligned) the contents of GMX_SIMD_REAL_WIDTH rvec elements sequentially into 3 SIMD registers
277 * The loaded output is:
278 * \p r0: { r[index][XX], r[index][YY], ... }
280 * \p r2: { ..., r[index+GMX_SIMD_REAL_WIDTH-1][YY], r[index+GMX_SIMD_REAL_WIDTH-1][ZZ] }
282 * \param[in] r Real to an rvec array, has to be aligned to SIMD register width
283 * \param[in] index Index of the first rvec triplet of reals to load
284 * \param[out] r0 Pointer to first SIMD register
285 * \param[out] r1 Pointer to second SIMD register
286 * \param[out] r2 Pointer to third SIMD register
288 static inline void simdLoadRvecs(const rvec *r,
294 const real *realPtr = r[index];
296 GMX_ASSERT(isSimdAligned(realPtr), "Pointer should be SIMD aligned");
298 *r0 = simdLoad(realPtr + 0*GMX_SIMD_REAL_WIDTH);
299 *r1 = simdLoad(realPtr + 1*GMX_SIMD_REAL_WIDTH);
300 *r2 = simdLoad(realPtr + 2*GMX_SIMD_REAL_WIDTH);
303 /*! \brief Store (aligned) 3 SIMD registers sequentially to GMX_SIMD_REAL_WIDTH rvec elements
305 * The stored output is:
306 * \p r[index] = { { r0[0], r0[1], ... }
308 * \p r[index+GMX_SIMD_REAL_WIDTH-1] = { ... , r2[GMX_SIMD_REAL_WIDTH-2], r2[GMX_SIMD_REAL_WIDTH-1] }
310 * \param[out] r Pointer to an rvec array, has to be aligned to SIMD register width
311 * \param[in] index Index of the first rvec triplet of reals to store to
312 * \param[in] r0 First SIMD register
313 * \param[in] r1 Second SIMD register
314 * \param[in] r2 Third SIMD register
316 static inline void simdStoreRvecs(rvec *r,
322 real *realPtr = r[index];
324 GMX_ASSERT(isSimdAligned(realPtr), "Pointer should be SIMD aligned");
326 store(realPtr + 0*GMX_SIMD_REAL_WIDTH, r0);
327 store(realPtr + 1*GMX_SIMD_REAL_WIDTH, r1);
328 store(realPtr + 2*GMX_SIMD_REAL_WIDTH, r2);
331 /*! \brief Integrate using leap-frog with single group T-scaling and SIMD
333 * \param[in] start Index of first atom to update
334 * \param[in] nrend Last atom to update: \p nrend - 1
335 * \param[in] dt The time step
336 * \param[in] invMass 1/mass per atom
337 * \param[in] tcstat Temperature coupling information
338 * \param[in] x Input coordinates
339 * \param[out] xprime Updated coordinates
340 * \param[inout] v Velocities
341 * \param[in] f Forces
344 updateMDLeapfrogSimpleSimd(int start,
347 const real * gmx_restrict invMass,
348 const t_grp_tcstat * tcstat,
349 const rvec * gmx_restrict x,
350 rvec * gmx_restrict xprime,
351 rvec * gmx_restrict v,
352 const rvec * gmx_restrict f)
354 SimdReal timestep(dt);
355 SimdReal lambdaSystem(tcstat[0].lambda);
357 /* We declare variables here, since code is often slower when declaring them inside the loop */
359 /* Note: We should implement a proper PaddedVector, so we don't need this check */
360 GMX_ASSERT(isSimdAligned(invMass), "invMass should be aligned");
362 for (int a = start; a < nrend; a += GMX_SIMD_REAL_WIDTH)
364 SimdReal invMass0, invMass1, invMass2;
365 expandScalarsToTriplets(simdLoad(invMass + a),
366 &invMass0, &invMass1, &invMass2);
370 simdLoadRvecs(v, a, &v0, &v1, &v2);
371 simdLoadRvecs(f, a, &f0, &f1, &f2);
373 v0 = fma(f0*invMass0, timestep, lambdaSystem*v0);
374 v1 = fma(f1*invMass1, timestep, lambdaSystem*v1);
375 v2 = fma(f2*invMass2, timestep, lambdaSystem*v2);
377 simdStoreRvecs(v, a, v0, v1, v2);
380 simdLoadRvecs(x, a, &x0, &x1, &x2);
382 SimdReal xprime0 = fma(v0, timestep, x0);
383 SimdReal xprime1 = fma(v1, timestep, x1);
384 SimdReal xprime2 = fma(v2, timestep, x2);
386 simdStoreRvecs(xprime, a, xprime0, xprime1, xprime2);
390 #endif // GMX_HAVE_SIMD_UPDATE
392 /*! \brief Sets the NEMD acceleration type */
393 enum class AccelerationType
398 /*! \brief Integrate using leap-frog with support for everything
400 * \tparam accelerationType Type of NEMD acceleration
401 * \param[in] start Index of first atom to update
402 * \param[in] nrend Last atom to update: \p nrend - 1
403 * \param[in] doNoseHoover If to apply Nose-Hoover T-coupling
404 * \param[in] dt The time step
405 * \param[in] dtPressureCouple Time step for pressure coupling, is 0 when no pressure coupling should be applied at this step
406 * \param[in] ir The input parameter record
407 * \param[in] md Atom properties
408 * \param[in] ekind Kinetic energy data
409 * \param[in] box The box dimensions
410 * \param[in] x Input coordinates
411 * \param[out] xprime Updated coordinates
412 * \param[inout] v Velocities
413 * \param[in] f Forces
414 * \param[in] nh_vxi Nose-Hoover velocity scaling factors
415 * \param[in] M Parrinello-Rahman scaling matrix
417 template<AccelerationType accelerationType>
419 updateMDLeapfrogGeneral(int start,
423 real dtPressureCouple,
424 const t_inputrec * ir,
425 const t_mdatoms * md,
426 const gmx_ekindata_t * ekind,
428 const rvec * gmx_restrict x,
429 rvec * gmx_restrict xprime,
430 rvec * gmx_restrict v,
431 const rvec * gmx_restrict f,
432 const double * gmx_restrict nh_vxi,
435 /* This is a version of the leap-frog integrator that supports
436 * all combinations of T-coupling, P-coupling and NEMD.
437 * Nose-Hoover uses the reversible leap-frog integrator from
438 * Holian et al. Phys Rev E 52(3) : 2338, 1995
441 const unsigned short * cTC = md->cTC;
442 const t_grp_tcstat * tcstat = ekind->tcstat;
444 const unsigned short * cACC = md->cACC;
445 const rvec * accel = ir->opts.acc;
446 const t_grp_acc * grpstat = ekind->grpstat;
448 const rvec * gmx_restrict invMassPerDim = md->invMassPerDim;
450 /* Initialize group values, changed later when multiple groups are used */
455 for (int n = start; n < nrend; n++)
461 real lg = tcstat[gt].lambda;
464 real cosineZ, vCosine;
465 #ifdef __INTEL_COMPILER
466 #pragma warning( disable : 280 )
468 switch (accelerationType)
470 case AccelerationType::none:
471 copy_rvec(v[n], vRel);
473 case AccelerationType::group:
478 /* Avoid scaling the group velocity */
479 rvec_sub(v[n], grpstat[ga].u, vRel);
481 case AccelerationType::cosine:
482 cosineZ = std::cos(x[n][ZZ]*static_cast<real>(M_PI)/box[ZZ][ZZ]);
483 vCosine = cosineZ*ekind->cosacc.vcos;
484 /* Avoid scaling the cosine profile velocity */
485 copy_rvec(v[n], vRel);
492 /* Here we account for multiple time stepping, by increasing
493 * the Nose-Hoover correction by nsttcouple
495 factorNH = 0.5*ir->nsttcouple*dt*nh_vxi[gt];
498 for (int d = 0; d < DIM; d++)
501 (lg*vRel[d] + (f[n][d]*invMassPerDim[n][d]*dt - factorNH*vRel[d]
502 - dtPressureCouple*iprod(M[d], vRel)))/(1 + factorNH);
503 switch (accelerationType)
505 case AccelerationType::none:
507 case AccelerationType::group:
508 /* Add back the mean velocity and apply acceleration */
509 vNew += grpstat[ga].u[d] + accel[ga][d]*dt;
511 case AccelerationType::cosine:
514 /* Add back the mean velocity and apply acceleration */
515 vNew += vCosine + cosineZ*ekind->cosacc.cos_accel*dt;
520 xprime[n][d] = x[n][d] + vNew*dt;
525 /*! \brief Handles the Leap-frog MD x and v integration */
526 static void do_update_md(int start,
530 const t_inputrec * ir,
531 const t_mdatoms * md,
532 const gmx_ekindata_t * ekind,
534 const rvec * gmx_restrict x,
535 rvec * gmx_restrict xprime,
536 rvec * gmx_restrict v,
537 const rvec * gmx_restrict f,
538 const double * gmx_restrict nh_vxi,
541 GMX_ASSERT(nrend == start || xprime != x, "For SIMD optimization certain compilers need to have xprime != x");
543 /* Note: Berendsen pressure scaling is handled after do_update_md() */
544 bool doTempCouple = isTemperatureCouplingStep(step, ir);
545 bool doNoseHoover = (ir->etc == etcNOSEHOOVER && doTempCouple);
546 bool doParrinelloRahman = (ir->epc == epcPARRINELLORAHMAN && isPressureCouplingStep(step, ir));
547 bool doPROffDiagonal = (doParrinelloRahman && (M[YY][XX] != 0 || M[ZZ][XX] != 0 || M[ZZ][YY] != 0));
549 real dtPressureCouple = (doParrinelloRahman ? ir->nstpcouple*dt : 0);
551 /* NEMD (also cosine) acceleration is applied in updateMDLeapFrogGeneral */
552 bool doAcceleration = (ekind->bNEMD || ekind->cosacc.cos_accel != 0);
554 if (doNoseHoover || doPROffDiagonal || doAcceleration)
557 if (!doParrinelloRahman)
559 /* We should not apply PR scaling at this step */
569 updateMDLeapfrogGeneral<AccelerationType::none>
570 (start, nrend, doNoseHoover, dt, dtPressureCouple,
571 ir, md, ekind, box, x, xprime, v, f, nh_vxi, stepM);
573 else if (ekind->bNEMD)
575 updateMDLeapfrogGeneral<AccelerationType::group>
576 (start, nrend, doNoseHoover, dt, dtPressureCouple,
577 ir, md, ekind, box, x, xprime, v, f, nh_vxi, stepM);
581 updateMDLeapfrogGeneral<AccelerationType::cosine>
582 (start, nrend, doNoseHoover, dt, dtPressureCouple,
583 ir, md, ekind, box, x, xprime, v, f, nh_vxi, stepM);
588 /* Use a simple and thus more efficient integration loop. */
589 /* The simple loop does not check for particle type (so it can
590 * be vectorized), which means we need to clear the velocities
591 * of virtual sites in advance, when present. Note that vsite
592 * velocities are computed after update and constraints from
593 * their displacement.
597 /* Note: The overhead of this loop is completely neligible */
598 clearVsiteVelocities(start, nrend, md->ptype, v);
601 /* We determine if we have a single T-coupling lambda value for all
602 * atoms. That allows for better SIMD acceleration in the template.
603 * If we do not do temperature coupling (in the run or this step),
604 * all scaling values are 1, so we effectively have a single value.
606 bool haveSingleTempScaleValue = (!doTempCouple || ekind->ngtc == 1);
608 /* Extract some pointers needed by all cases */
609 const unsigned short *cTC = md->cTC;
610 const t_grp_tcstat *tcstat = ekind->tcstat;
611 const rvec *invMassPerDim = md->invMassPerDim;
613 if (doParrinelloRahman)
615 GMX_ASSERT(!doPROffDiagonal, "updateMDLeapfrogSimple only support diagonal Parrinello-Rahman scaling matrices");
618 for (int d = 0; d < DIM; d++)
623 if (haveSingleTempScaleValue)
625 updateMDLeapfrogSimple
626 <NumTempScaleValues::single,
627 ApplyParrinelloRahmanVScaling::diagonal>
628 (start, nrend, dt, dtPressureCouple,
629 invMassPerDim, tcstat, cTC, diagM, x, xprime, v, f);
633 updateMDLeapfrogSimple
634 <NumTempScaleValues::multiple,
635 ApplyParrinelloRahmanVScaling::diagonal>
636 (start, nrend, dt, dtPressureCouple,
637 invMassPerDim, tcstat, cTC, diagM, x, xprime, v, f);
642 if (haveSingleTempScaleValue)
644 /* Note that modern compilers are pretty good at vectorizing
645 * updateMDLeapfrogSimple(). But the SIMD version will still
646 * be faster because invMass lowers the cache pressure
647 * compared to invMassPerDim.
649 #if GMX_HAVE_SIMD_UPDATE
650 /* Check if we can use invmass instead of invMassPerDim */
651 if (!md->havePartiallyFrozenAtoms)
653 updateMDLeapfrogSimpleSimd
654 (start, nrend, dt, md->invmass, tcstat, x, xprime, v, f);
659 updateMDLeapfrogSimple
660 <NumTempScaleValues::single,
661 ApplyParrinelloRahmanVScaling::no>
662 (start, nrend, dt, dtPressureCouple,
663 invMassPerDim, tcstat, cTC, nullptr, x, xprime, v, f);
668 updateMDLeapfrogSimple
669 <NumTempScaleValues::multiple,
670 ApplyParrinelloRahmanVScaling::no>
671 (start, nrend, dt, dtPressureCouple,
672 invMassPerDim, tcstat, cTC, nullptr, x, xprime, v, f);
678 static void do_update_vv_vel(int start, int nrend, real dt,
679 const rvec accel[], const ivec nFreeze[], const real invmass[],
680 const unsigned short ptype[], const unsigned short cFREEZE[],
681 const unsigned short cACC[], rvec v[], const rvec f[],
682 gmx_bool bExtended, real veta, real alpha)
690 g = 0.25*dt*veta*alpha;
692 mv2 = gmx::series_sinhx(g);
699 for (n = start; n < nrend; n++)
701 real w_dt = invmass[n]*dt;
711 for (d = 0; d < DIM; d++)
713 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
715 v[n][d] = mv1*(mv1*v[n][d] + 0.5*(w_dt*mv2*f[n][d]))+0.5*accel[ga][d]*dt;
723 } /* do_update_vv_vel */
725 static void do_update_vv_pos(int start, int nrend, real dt,
726 const ivec nFreeze[],
727 const unsigned short ptype[], const unsigned short cFREEZE[],
728 const rvec x[], rvec xprime[], const rvec v[],
729 gmx_bool bExtended, real veta)
735 /* Would it make more sense if Parrinello-Rahman was put here? */
740 mr2 = gmx::series_sinhx(g);
748 for (n = start; n < nrend; n++)
756 for (d = 0; d < DIM; d++)
758 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
760 xprime[n][d] = mr1*(mr1*x[n][d]+mr2*dt*v[n][d]);
764 xprime[n][d] = x[n][d];
768 } /* do_update_vv_pos */
770 gmx_stochd_t::gmx_stochd_t(const t_inputrec *ir)
772 const t_grpopts *opts = &ir->opts;
773 const int ngtc = opts->ngtc;
779 else if (EI_SD(ir->eI))
784 for (int gt = 0; gt < ngtc; gt++)
786 if (opts->tau_t[gt] > 0)
788 sdc[gt].em = std::exp(-ir->delta_t/opts->tau_t[gt]);
792 /* No friction and noise on this group */
797 else if (ETC_ANDERSEN(ir->etc))
799 randomize_group.resize(ngtc);
800 boltzfac.resize(ngtc);
802 /* for now, assume that all groups, if randomized, are randomized at the same rate, i.e. tau_t is the same. */
803 /* since constraint groups don't necessarily match up with temperature groups! This is checked in readir.c */
805 for (int gt = 0; gt < ngtc; gt++)
807 real reft = std::max<real>(0, opts->ref_t[gt]);
808 if ((opts->tau_t[gt] > 0) && (reft > 0)) /* tau_t or ref_t = 0 means that no randomization is done */
810 randomize_group[gt] = true;
811 boltzfac[gt] = BOLTZ*opts->ref_t[gt];
815 randomize_group[gt] = false;
821 void update_temperature_constants(gmx_update_t *upd, const t_inputrec *ir)
825 if (ir->bd_fric != 0)
827 for (int gt = 0; gt < ir->opts.ngtc; gt++)
829 upd->sd->bd_rf[gt] = std::sqrt(2.0*BOLTZ*ir->opts.ref_t[gt]/(ir->bd_fric*ir->delta_t));
834 for (int gt = 0; gt < ir->opts.ngtc; gt++)
836 upd->sd->bd_rf[gt] = std::sqrt(2.0*BOLTZ*ir->opts.ref_t[gt]);
842 for (int gt = 0; gt < ir->opts.ngtc; gt++)
844 real kT = BOLTZ*ir->opts.ref_t[gt];
845 /* The mass is accounted for later, since this differs per atom */
846 upd->sd->sdsig[gt].V = std::sqrt(kT*(1 - upd->sd->sdc[gt].em*upd->sd->sdc[gt].em));
851 gmx_update_t *init_update(const t_inputrec *ir,
852 gmx::BoxDeformation *deform)
854 gmx_update_t *upd = new(gmx_update_t);
856 upd->sd = gmx::compat::make_unique<gmx_stochd_t>(ir);
858 update_temperature_constants(upd, ir);
860 upd->xp.resizeWithPadding(0);
862 upd->deform = deform;
867 void update_realloc(gmx_update_t *upd, int natoms)
869 GMX_ASSERT(upd, "upd must be allocated before its fields can be reallocated");
871 upd->xp.resizeWithPadding(natoms);
874 /*! \brief Sets the SD update type */
875 enum class SDUpdate : int
877 ForcesOnly, FrictionAndNoiseOnly, Combined
880 /*! \brief SD integrator update
882 * Two phases are required in the general case of a constrained
883 * update, the first phase from the contribution of forces, before
884 * applying constraints, and then a second phase applying the friction
885 * and noise, and then further constraining. For details, see
888 * Without constraints, the two phases can be combined, for
891 * Thus three instantiations of this templated function will be made,
892 * two with only one contribution, and one with both contributions. */
893 template <SDUpdate updateType>
895 doSDUpdateGeneral(const gmx_stochd_t &sd,
896 int start, int nrend, real dt,
897 const rvec accel[], const ivec nFreeze[],
898 const real invmass[], const unsigned short ptype[],
899 const unsigned short cFREEZE[], const unsigned short cACC[],
900 const unsigned short cTC[],
901 const rvec x[], rvec xprime[], rvec v[], const rvec f[],
902 int64_t step, int seed, const int *gatindex)
904 // cTC, cACC and cFreeze can be nullptr any time, but various
905 // instantiations do not make sense with particular pointer
907 if (updateType == SDUpdate::ForcesOnly)
909 GMX_ASSERT(f != nullptr, "SD update with only forces requires forces");
910 GMX_ASSERT(cTC == nullptr, "SD update with only forces cannot handle temperature groups");
912 if (updateType == SDUpdate::FrictionAndNoiseOnly)
914 GMX_ASSERT(f == nullptr, "SD update with only noise cannot handle forces");
915 GMX_ASSERT(cACC == nullptr, "SD update with only noise cannot handle acceleration groups");
917 if (updateType == SDUpdate::Combined)
919 GMX_ASSERT(f != nullptr, "SD update with forces and noise requires forces");
922 // Even 0 bits internal counter gives 2x64 ints (more than enough for three table lookups)
923 gmx::ThreeFry2x64<0> rng(seed, gmx::RandomDomain::UpdateCoordinates);
924 gmx::TabulatedNormalDistribution<real, 14> dist;
926 for (int n = start; n < nrend; n++)
928 int globalAtomIndex = gatindex ? gatindex[n] : n;
929 rng.restart(step, globalAtomIndex);
932 real inverseMass = invmass[n];
933 real invsqrtMass = std::sqrt(inverseMass);
935 int freezeGroup = cFREEZE ? cFREEZE[n] : 0;
936 int accelerationGroup = cACC ? cACC[n] : 0;
937 int temperatureGroup = cTC ? cTC[n] : 0;
939 for (int d = 0; d < DIM; d++)
941 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[freezeGroup][d])
943 if (updateType == SDUpdate::ForcesOnly)
945 real vn = v[n][d] + (inverseMass*f[n][d] + accel[accelerationGroup][d])*dt;
947 // Simple position update.
948 xprime[n][d] = x[n][d] + v[n][d]*dt;
950 else if (updateType == SDUpdate::FrictionAndNoiseOnly)
953 v[n][d] = (vn*sd.sdc[temperatureGroup].em +
954 invsqrtMass*sd.sdsig[temperatureGroup].V*dist(rng));
955 // The previous phase already updated the
956 // positions with a full v*dt term that must
957 // now be half removed.
958 xprime[n][d] = xprime[n][d] + 0.5*(v[n][d] - vn)*dt;
962 real vn = v[n][d] + (inverseMass*f[n][d] + accel[accelerationGroup][d])*dt;
963 v[n][d] = (vn*sd.sdc[temperatureGroup].em +
964 invsqrtMass*sd.sdsig[temperatureGroup].V*dist(rng));
965 // Here we include half of the friction+noise
966 // update of v into the position update.
967 xprime[n][d] = x[n][d] + 0.5*(vn + v[n][d])*dt;
972 // When using constraints, the update is split into
973 // two phases, but we only need to zero the update of
974 // virtual, shell or frozen particles in at most one
976 if (updateType != SDUpdate::FrictionAndNoiseOnly)
979 xprime[n][d] = x[n][d];
986 static void do_update_bd(int start, int nrend, real dt,
987 const ivec nFreeze[],
988 const real invmass[], const unsigned short ptype[],
989 const unsigned short cFREEZE[], const unsigned short cTC[],
990 const rvec x[], rvec xprime[], rvec v[],
991 const rvec f[], real friction_coefficient,
992 const real *rf, int64_t step, int seed,
995 /* note -- these appear to be full step velocities . . . */
1000 // Use 1 bit of internal counters to give us 2*2 64-bits values per stream
1001 // Each 64-bit value is enough for 4 normal distribution table numbers.
1002 gmx::ThreeFry2x64<0> rng(seed, gmx::RandomDomain::UpdateCoordinates);
1003 gmx::TabulatedNormalDistribution<real, 14> dist;
1005 if (friction_coefficient != 0)
1007 invfr = 1.0/friction_coefficient;
1010 for (n = start; (n < nrend); n++)
1012 int ng = gatindex ? gatindex[n] : n;
1014 rng.restart(step, ng);
1025 for (d = 0; (d < DIM); d++)
1027 if ((ptype[n] != eptVSite) && (ptype[n] != eptShell) && !nFreeze[gf][d])
1029 if (friction_coefficient != 0)
1031 vn = invfr*f[n][d] + rf[gt]*dist(rng);
1035 /* NOTE: invmass = 2/(mass*friction_constant*dt) */
1036 vn = 0.5*invmass[n]*f[n][d]*dt
1037 + std::sqrt(0.5*invmass[n])*rf[gt]*dist(rng);
1041 xprime[n][d] = x[n][d]+vn*dt;
1046 xprime[n][d] = x[n][d];
1052 static void calc_ke_part_normal(const rvec v[], const t_grpopts *opts, const t_mdatoms *md,
1053 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel)
1056 t_grp_tcstat *tcstat = ekind->tcstat;
1057 t_grp_acc *grpstat = ekind->grpstat;
1058 int nthread, thread;
1060 /* three main: VV with AveVel, vv with AveEkin, leap with AveEkin. Leap with AveVel is also
1061 an option, but not supported now.
1062 bEkinAveVel: If TRUE, we sum into ekin, if FALSE, into ekinh.
1065 /* group velocities are calculated in update_ekindata and
1066 * accumulated in acumulate_groups.
1067 * Now the partial global and groups ekin.
1069 for (g = 0; (g < opts->ngtc); g++)
1071 copy_mat(tcstat[g].ekinh, tcstat[g].ekinh_old);
1074 clear_mat(tcstat[g].ekinf);
1075 tcstat[g].ekinscalef_nhc = 1.0; /* need to clear this -- logic is complicated! */
1079 clear_mat(tcstat[g].ekinh);
1082 ekind->dekindl_old = ekind->dekindl;
1083 nthread = gmx_omp_nthreads_get(emntUpdate);
1085 #pragma omp parallel for num_threads(nthread) schedule(static)
1086 for (thread = 0; thread < nthread; thread++)
1088 // This OpenMP only loops over arrays and does not call any functions
1089 // or memory allocation. It should not be able to throw, so for now
1090 // we do not need a try/catch wrapper.
1091 int start_t, end_t, n;
1099 start_t = ((thread+0)*md->homenr)/nthread;
1100 end_t = ((thread+1)*md->homenr)/nthread;
1102 ekin_sum = ekind->ekin_work[thread];
1103 dekindl_sum = ekind->dekindl_work[thread];
1105 for (gt = 0; gt < opts->ngtc; gt++)
1107 clear_mat(ekin_sum[gt]);
1113 for (n = start_t; n < end_t; n++)
1123 hm = 0.5*md->massT[n];
1125 for (d = 0; (d < DIM); d++)
1127 v_corrt[d] = v[n][d] - grpstat[ga].u[d];
1129 for (d = 0; (d < DIM); d++)
1131 for (m = 0; (m < DIM); m++)
1133 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1134 ekin_sum[gt][m][d] += hm*v_corrt[m]*v_corrt[d];
1137 if (md->nMassPerturbed && md->bPerturbed[n])
1140 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1146 for (thread = 0; thread < nthread; thread++)
1148 for (g = 0; g < opts->ngtc; g++)
1152 m_add(tcstat[g].ekinf, ekind->ekin_work[thread][g],
1157 m_add(tcstat[g].ekinh, ekind->ekin_work[thread][g],
1162 ekind->dekindl += *ekind->dekindl_work[thread];
1165 inc_nrnb(nrnb, eNR_EKIN, md->homenr);
1168 static void calc_ke_part_visc(const matrix box, const rvec x[], const rvec v[],
1169 const t_grpopts *opts, const t_mdatoms *md,
1170 gmx_ekindata_t *ekind,
1171 t_nrnb *nrnb, gmx_bool bEkinAveVel)
1173 int start = 0, homenr = md->homenr;
1174 int g, d, n, m, gt = 0;
1177 t_grp_tcstat *tcstat = ekind->tcstat;
1178 t_cos_acc *cosacc = &(ekind->cosacc);
1183 for (g = 0; g < opts->ngtc; g++)
1185 copy_mat(ekind->tcstat[g].ekinh, ekind->tcstat[g].ekinh_old);
1186 clear_mat(ekind->tcstat[g].ekinh);
1188 ekind->dekindl_old = ekind->dekindl;
1190 fac = 2*M_PI/box[ZZ][ZZ];
1193 for (n = start; n < start+homenr; n++)
1199 hm = 0.5*md->massT[n];
1201 /* Note that the times of x and v differ by half a step */
1202 /* MRS -- would have to be changed for VV */
1203 cosz = std::cos(fac*x[n][ZZ]);
1204 /* Calculate the amplitude of the new velocity profile */
1205 mvcos += 2*cosz*md->massT[n]*v[n][XX];
1207 copy_rvec(v[n], v_corrt);
1208 /* Subtract the profile for the kinetic energy */
1209 v_corrt[XX] -= cosz*cosacc->vcos;
1210 for (d = 0; (d < DIM); d++)
1212 for (m = 0; (m < DIM); m++)
1214 /* if we're computing a full step velocity, v_corrt[d] has v(t). Otherwise, v(t+dt/2) */
1217 tcstat[gt].ekinf[m][d] += hm*v_corrt[m]*v_corrt[d];
1221 tcstat[gt].ekinh[m][d] += hm*v_corrt[m]*v_corrt[d];
1225 if (md->nPerturbed && md->bPerturbed[n])
1227 /* The minus sign here might be confusing.
1228 * The kinetic contribution from dH/dl doesn't come from
1229 * d m(l)/2 v^2 / dl, but rather from d p^2/2m(l) / dl,
1230 * where p are the momenta. The difference is only a minus sign.
1232 dekindl -= 0.5*(md->massB[n] - md->massA[n])*iprod(v_corrt, v_corrt);
1235 ekind->dekindl = dekindl;
1236 cosacc->mvcos = mvcos;
1238 inc_nrnb(nrnb, eNR_EKIN, homenr);
1241 void calc_ke_part(const t_state *state, const t_grpopts *opts, const t_mdatoms *md,
1242 gmx_ekindata_t *ekind, t_nrnb *nrnb, gmx_bool bEkinAveVel)
1244 if (ekind->cosacc.cos_accel == 0)
1246 calc_ke_part_normal(state->v.rvec_array(), opts, md, ekind, nrnb, bEkinAveVel);
1250 calc_ke_part_visc(state->box, state->x.rvec_array(), state->v.rvec_array(), opts, md, ekind, nrnb, bEkinAveVel);
1254 extern void init_ekinstate(ekinstate_t *ekinstate, const t_inputrec *ir)
1256 ekinstate->ekin_n = ir->opts.ngtc;
1257 snew(ekinstate->ekinh, ekinstate->ekin_n);
1258 snew(ekinstate->ekinf, ekinstate->ekin_n);
1259 snew(ekinstate->ekinh_old, ekinstate->ekin_n);
1260 ekinstate->ekinscalef_nhc.resize(ekinstate->ekin_n);
1261 ekinstate->ekinscaleh_nhc.resize(ekinstate->ekin_n);
1262 ekinstate->vscale_nhc.resize(ekinstate->ekin_n);
1263 ekinstate->dekindl = 0;
1264 ekinstate->mvcos = 0;
1267 void update_ekinstate(ekinstate_t *ekinstate, const gmx_ekindata_t *ekind)
1271 for (i = 0; i < ekinstate->ekin_n; i++)
1273 copy_mat(ekind->tcstat[i].ekinh, ekinstate->ekinh[i]);
1274 copy_mat(ekind->tcstat[i].ekinf, ekinstate->ekinf[i]);
1275 copy_mat(ekind->tcstat[i].ekinh_old, ekinstate->ekinh_old[i]);
1276 ekinstate->ekinscalef_nhc[i] = ekind->tcstat[i].ekinscalef_nhc;
1277 ekinstate->ekinscaleh_nhc[i] = ekind->tcstat[i].ekinscaleh_nhc;
1278 ekinstate->vscale_nhc[i] = ekind->tcstat[i].vscale_nhc;
1281 copy_mat(ekind->ekin, ekinstate->ekin_total);
1282 ekinstate->dekindl = ekind->dekindl;
1283 ekinstate->mvcos = ekind->cosacc.mvcos;
1287 void restore_ekinstate_from_state(const t_commrec *cr,
1288 gmx_ekindata_t *ekind, const ekinstate_t *ekinstate)
1294 for (i = 0; i < ekinstate->ekin_n; i++)
1296 copy_mat(ekinstate->ekinh[i], ekind->tcstat[i].ekinh);
1297 copy_mat(ekinstate->ekinf[i], ekind->tcstat[i].ekinf);
1298 copy_mat(ekinstate->ekinh_old[i], ekind->tcstat[i].ekinh_old);
1299 ekind->tcstat[i].ekinscalef_nhc = ekinstate->ekinscalef_nhc[i];
1300 ekind->tcstat[i].ekinscaleh_nhc = ekinstate->ekinscaleh_nhc[i];
1301 ekind->tcstat[i].vscale_nhc = ekinstate->vscale_nhc[i];
1304 copy_mat(ekinstate->ekin_total, ekind->ekin);
1306 ekind->dekindl = ekinstate->dekindl;
1307 ekind->cosacc.mvcos = ekinstate->mvcos;
1308 n = ekinstate->ekin_n;
1313 gmx_bcast(sizeof(n), &n, cr);
1314 for (i = 0; i < n; i++)
1316 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh[0][0]),
1317 ekind->tcstat[i].ekinh[0], cr);
1318 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinf[0][0]),
1319 ekind->tcstat[i].ekinf[0], cr);
1320 gmx_bcast(DIM*DIM*sizeof(ekind->tcstat[i].ekinh_old[0][0]),
1321 ekind->tcstat[i].ekinh_old[0], cr);
1323 gmx_bcast(sizeof(ekind->tcstat[i].ekinscalef_nhc),
1324 &(ekind->tcstat[i].ekinscalef_nhc), cr);
1325 gmx_bcast(sizeof(ekind->tcstat[i].ekinscaleh_nhc),
1326 &(ekind->tcstat[i].ekinscaleh_nhc), cr);
1327 gmx_bcast(sizeof(ekind->tcstat[i].vscale_nhc),
1328 &(ekind->tcstat[i].vscale_nhc), cr);
1330 gmx_bcast(DIM*DIM*sizeof(ekind->ekin[0][0]),
1331 ekind->ekin[0], cr);
1333 gmx_bcast(sizeof(ekind->dekindl), &ekind->dekindl, cr);
1334 gmx_bcast(sizeof(ekind->cosacc.mvcos), &ekind->cosacc.mvcos, cr);
1338 void update_tcouple(int64_t step,
1339 const t_inputrec *inputrec,
1341 gmx_ekindata_t *ekind,
1342 const t_extmass *MassQ,
1343 const t_mdatoms *md)
1346 bool doTemperatureCoupling = false;
1348 /* if using vv with trotter decomposition methods, we do this elsewhere in the code */
1349 if (!(EI_VV(inputrec->eI) &&
1350 (inputrecNvtTrotter(inputrec) || inputrecNptTrotter(inputrec) || inputrecNphTrotter(inputrec))))
1352 doTemperatureCoupling = isTemperatureCouplingStep(step, inputrec);
1355 if (doTemperatureCoupling)
1357 real dttc = inputrec->nsttcouple*inputrec->delta_t;
1359 switch (inputrec->etc)
1364 berendsen_tcoupl(inputrec, ekind, dttc, state->therm_integral);
1367 nosehoover_tcoupl(&(inputrec->opts), ekind, dttc,
1368 state->nosehoover_xi.data(), state->nosehoover_vxi.data(), MassQ);
1371 vrescale_tcoupl(inputrec, step, ekind, dttc,
1372 state->therm_integral.data());
1375 /* rescale in place here */
1376 if (EI_VV(inputrec->eI))
1378 rescale_velocities(ekind, md, 0, md->homenr, state->v.rvec_array());
1383 /* Set the T scaling lambda to 1 to have no scaling */
1384 for (int i = 0; (i < inputrec->opts.ngtc); i++)
1386 ekind->tcstat[i].lambda = 1.0;
1391 /*! \brief Computes the atom range for a thread to operate on, ensured SIMD aligned ranges
1393 * \param[in] numThreads The number of threads to divide atoms over
1394 * \param[in] threadIndex The thread to get the range for
1395 * \param[in] numAtoms The total number of atoms (on this rank)
1396 * \param[out] startAtom The start of the atom range
1397 * \param[out] endAtom The end of the atom range, note that this is in general not a multiple of the SIMD width
1399 static void getThreadAtomRange(int numThreads, int threadIndex, int numAtoms,
1400 int *startAtom, int *endAtom)
1402 #if GMX_HAVE_SIMD_UPDATE
1403 constexpr int blockSize = GMX_SIMD_REAL_WIDTH;
1405 constexpr int blockSize = 1;
1407 int numBlocks = (numAtoms + blockSize - 1)/blockSize;
1409 *startAtom = ((numBlocks* threadIndex )/numThreads)*blockSize;
1410 *endAtom = ((numBlocks*(threadIndex + 1))/numThreads)*blockSize;
1411 if (threadIndex == numThreads - 1)
1413 *endAtom = numAtoms;
1417 void update_pcouple_before_coordinates(FILE *fplog,
1419 const t_inputrec *inputrec,
1421 matrix parrinellorahmanMu,
1425 /* Berendsen P-coupling is completely handled after the coordinate update.
1426 * Trotter P-coupling is handled by separate calls to trotter_update().
1428 if (inputrec->epc == epcPARRINELLORAHMAN &&
1429 isPressureCouplingStep(step, inputrec))
1431 real dtpc = inputrec->nstpcouple*inputrec->delta_t;
1433 parrinellorahman_pcoupl(fplog, step, inputrec, dtpc, state->pres_prev,
1434 state->box, state->box_rel, state->boxv,
1435 M, parrinellorahmanMu, bInitStep);
1439 void constrain_velocities(int64_t step,
1440 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1443 gmx::Constraints *constr,
1455 * APPLY CONSTRAINTS:
1462 /* clear out constraints before applying */
1463 clear_mat(vir_part);
1465 /* Constrain the coordinates upd->xp */
1466 constr->apply(do_log, do_ene,
1468 state->x.rvec_array(), state->v.rvec_array(), state->v.rvec_array(),
1470 state->lambda[efptBONDED], dvdlambda,
1471 nullptr, bCalcVir ? &vir_con : nullptr, ConstraintVariable::Velocities);
1475 m_add(vir_part, vir_con, vir_part);
1480 void constrain_coordinates(int64_t step,
1481 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1485 gmx::Constraints *constr,
1498 /* clear out constraints before applying */
1499 clear_mat(vir_part);
1501 /* Constrain the coordinates upd->xp */
1502 constr->apply(do_log, do_ene,
1504 state->x.rvec_array(), upd->xp.rvec_array(), nullptr,
1506 state->lambda[efptBONDED], dvdlambda,
1507 as_rvec_array(state->v.data()), bCalcVir ? &vir_con : nullptr, ConstraintVariable::Positions);
1511 m_add(vir_part, vir_con, vir_part);
1517 update_sd_second_half(int64_t step,
1518 real *dvdlambda, /* the contribution to be added to the bonded interactions */
1519 const t_inputrec *inputrec, /* input record and box stuff */
1520 const t_mdatoms *md,
1522 const t_commrec *cr,
1524 gmx_wallcycle_t wcycle,
1526 gmx::Constraints *constr,
1534 if (inputrec->eI == eiSD1)
1537 int homenr = md->homenr;
1539 /* Cast delta_t from double to real to make the integrators faster.
1540 * The only reason for having delta_t double is to get accurate values
1541 * for t=delta_t*step when step is larger than float precision.
1542 * For integration dt the accuracy of real suffices, since with
1543 * integral += dt*integrand the increment is nearly always (much) smaller
1544 * than the integral (and the integrand has real precision).
1546 real dt = inputrec->delta_t;
1548 wallcycle_start(wcycle, ewcUPDATE);
1550 nth = gmx_omp_nthreads_get(emntUpdate);
1552 #pragma omp parallel for num_threads(nth) schedule(static)
1553 for (th = 0; th < nth; th++)
1557 int start_th, end_th;
1558 getThreadAtomRange(nth, th, homenr, &start_th, &end_th);
1560 doSDUpdateGeneral<SDUpdate::FrictionAndNoiseOnly>
1562 start_th, end_th, dt,
1563 inputrec->opts.acc, inputrec->opts.nFreeze,
1564 md->invmass, md->ptype,
1565 md->cFREEZE, nullptr, md->cTC,
1566 state->x.rvec_array(), upd->xp.rvec_array(),
1567 state->v.rvec_array(), nullptr,
1568 step, inputrec->ld_seed,
1569 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1571 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1573 inc_nrnb(nrnb, eNR_UPDATE, homenr);
1574 wallcycle_stop(wcycle, ewcUPDATE);
1576 /* Constrain the coordinates upd->xp for half a time step */
1577 constr->apply(do_log, do_ene,
1579 state->x.rvec_array(), upd->xp.rvec_array(), nullptr,
1581 state->lambda[efptBONDED], dvdlambda,
1582 as_rvec_array(state->v.data()), nullptr, ConstraintVariable::Positions);
1586 void finish_update(const t_inputrec *inputrec, /* input record and box stuff */
1587 const t_mdatoms *md,
1589 const t_graph *graph,
1591 gmx_wallcycle_t wcycle,
1593 const gmx::Constraints *constr)
1595 int homenr = md->homenr;
1597 /* We must always unshift after updating coordinates; if we did not shake
1598 x was shifted in do_force */
1600 /* NOTE Currently we always integrate to a temporary buffer and
1601 * then copy the results back. */
1603 wallcycle_start_nocount(wcycle, ewcUPDATE);
1605 if (md->cFREEZE != nullptr && constr != nullptr)
1607 /* If we have atoms that are frozen along some, but not all
1608 * dimensions, then any constraints will have moved them also along
1609 * the frozen dimensions. To freeze such degrees of freedom
1610 * we copy them back here to later copy them forward. It would
1611 * be more elegant and slightly more efficient to copies zero
1612 * times instead of twice, but the graph case below prevents this.
1614 const ivec *nFreeze = inputrec->opts.nFreeze;
1615 bool partialFreezeAndConstraints = false;
1616 for (int g = 0; g < inputrec->opts.ngfrz; g++)
1618 int numFreezeDim = nFreeze[g][XX] + nFreeze[g][YY] + nFreeze[g][ZZ];
1619 if (numFreezeDim > 0 && numFreezeDim < 3)
1621 partialFreezeAndConstraints = true;
1624 if (partialFreezeAndConstraints)
1626 auto xp = makeArrayRef(upd->xp).subArray(0, homenr);
1627 auto x = makeConstArrayRef(state->x).subArray(0, homenr);
1628 for (int i = 0; i < homenr; i++)
1630 int g = md->cFREEZE[i];
1632 for (int d = 0; d < DIM; d++)
1643 if (graph && (graph->nnodes > 0))
1645 unshift_x(graph, state->box, state->x.rvec_array(), upd->xp.rvec_array());
1646 if (TRICLINIC(state->box))
1648 inc_nrnb(nrnb, eNR_SHIFTX, 2*graph->nnodes);
1652 inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
1657 auto xp = makeConstArrayRef(upd->xp).subArray(0, homenr);
1658 auto x = makeArrayRef(state->x).subArray(0, homenr);
1659 #ifndef __clang_analyzer__
1660 int gmx_unused nth = gmx_omp_nthreads_get(emntUpdate);
1662 #pragma omp parallel for num_threads(nth) schedule(static)
1663 for (int i = 0; i < homenr; i++)
1665 // Trivial statement, does not throw
1669 wallcycle_stop(wcycle, ewcUPDATE);
1671 /* ############# END the update of velocities and positions ######### */
1674 void update_pcouple_after_coordinates(FILE *fplog,
1676 const t_inputrec *inputrec,
1677 const t_mdatoms *md,
1678 const matrix pressure,
1679 const matrix forceVirial,
1680 const matrix constraintVirial,
1681 const matrix parrinellorahmanMu,
1687 int homenr = md->homenr;
1689 /* Cast to real for faster code, no loss in precision (see comment above) */
1690 real dt = inputrec->delta_t;
1693 /* now update boxes */
1694 switch (inputrec->epc)
1698 case (epcBERENDSEN):
1699 if (isPressureCouplingStep(step, inputrec))
1701 real dtpc = inputrec->nstpcouple*dt;
1703 berendsen_pcoupl(fplog, step, inputrec, dtpc,
1704 pressure, state->box,
1705 forceVirial, constraintVirial,
1706 mu, &state->baros_integral);
1707 berendsen_pscale(inputrec, mu, state->box, state->box_rel,
1708 start, homenr, state->x.rvec_array(),
1712 case (epcPARRINELLORAHMAN):
1713 if (isPressureCouplingStep(step, inputrec))
1715 /* The box velocities were updated in do_pr_pcoupl,
1716 * but we dont change the box vectors until we get here
1717 * since we need to be able to shift/unshift above.
1719 real dtpc = inputrec->nstpcouple*dt;
1720 for (int i = 0; i < DIM; i++)
1722 for (int m = 0; m <= i; m++)
1724 state->box[i][m] += dtpc*state->boxv[i][m];
1727 preserve_box_shape(inputrec, state->box_rel, state->box);
1729 /* Scale the coordinates */
1730 auto x = state->x.rvec_array();
1731 for (int n = start; n < start + homenr; n++)
1733 tmvmul_ur0(parrinellorahmanMu, x[n], x[n]);
1738 switch (inputrec->epct)
1740 case (epctISOTROPIC):
1741 /* DIM * eta = ln V. so DIM*eta_new = DIM*eta_old + DIM*dt*veta =>
1742 ln V_new = ln V_old + 3*dt*veta => V_new = V_old*exp(3*dt*veta) =>
1743 Side length scales as exp(veta*dt) */
1745 msmul(state->box, std::exp(state->veta*dt), state->box);
1747 /* Relate veta to boxv. veta = d(eta)/dT = (1/DIM)*1/V dV/dT.
1748 o If we assume isotropic scaling, and box length scaling
1749 factor L, then V = L^DIM (det(M)). So dV/dt = DIM
1750 L^(DIM-1) dL/dt det(M), and veta = (1/L) dL/dt. The
1751 determinant of B is L^DIM det(M), and the determinant
1752 of dB/dt is (dL/dT)^DIM det (M). veta will be
1753 (det(dB/dT)/det(B))^(1/3). Then since M =
1754 B_new*(vol_new)^(1/3), dB/dT_new = (veta_new)*B(new). */
1756 msmul(state->box, state->veta, state->boxv);
1768 auto localX = makeArrayRef(state->x).subArray(start, homenr);
1769 upd->deform->apply(localX, state->box, step);
1773 void update_coords(int64_t step,
1774 const t_inputrec *inputrec, /* input record and box stuff */
1775 const t_mdatoms *md,
1777 gmx::ArrayRefWithPadding<gmx::RVec> f,
1778 const t_fcdata *fcd,
1779 const gmx_ekindata_t *ekind,
1783 const t_commrec *cr, /* these shouldn't be here -- need to think about it */
1784 const gmx::Constraints *constr)
1786 gmx_bool bDoConstr = (nullptr != constr);
1788 /* Running the velocity half does nothing except for velocity verlet */
1789 if ((UpdatePart == etrtVELOCITY1 || UpdatePart == etrtVELOCITY2) &&
1790 !EI_VV(inputrec->eI))
1792 gmx_incons("update_coords called for velocity without VV integrator");
1795 int homenr = md->homenr;
1797 /* Cast to real for faster code, no loss in precision (see comment above) */
1798 real dt = inputrec->delta_t;
1800 /* We need to update the NMR restraint history when time averaging is used */
1801 if (state->flags & (1<<estDISRE_RM3TAV))
1803 update_disres_history(fcd, &state->hist);
1805 if (state->flags & (1<<estORIRE_DTAV))
1807 update_orires_history(fcd, &state->hist);
1810 /* ############# START The update of velocities and positions ######### */
1811 int nth = gmx_omp_nthreads_get(emntUpdate);
1813 #pragma omp parallel for num_threads(nth) schedule(static)
1814 for (int th = 0; th < nth; th++)
1818 int start_th, end_th;
1819 getThreadAtomRange(nth, th, homenr, &start_th, &end_th);
1821 const rvec *x_rvec = state->x.rvec_array();
1822 rvec *xp_rvec = upd->xp.rvec_array();
1823 rvec *v_rvec = state->v.rvec_array();
1824 const rvec *f_rvec = as_rvec_array(f.unpaddedArrayRef().data());
1826 switch (inputrec->eI)
1829 do_update_md(start_th, end_th, step, dt,
1830 inputrec, md, ekind, state->box,
1831 x_rvec, xp_rvec, v_rvec, f_rvec,
1832 state->nosehoover_vxi.data(), M);
1837 // With constraints, the SD update is done in 2 parts
1838 doSDUpdateGeneral<SDUpdate::ForcesOnly>
1840 start_th, end_th, dt,
1841 inputrec->opts.acc, inputrec->opts.nFreeze,
1842 md->invmass, md->ptype,
1843 md->cFREEZE, md->cACC, nullptr,
1844 x_rvec, xp_rvec, v_rvec, f_rvec,
1845 step, inputrec->ld_seed, nullptr);
1849 doSDUpdateGeneral<SDUpdate::Combined>
1851 start_th, end_th, dt,
1852 inputrec->opts.acc, inputrec->opts.nFreeze,
1853 md->invmass, md->ptype,
1854 md->cFREEZE, md->cACC, md->cTC,
1855 x_rvec, xp_rvec, v_rvec, f_rvec,
1856 step, inputrec->ld_seed,
1857 DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1861 do_update_bd(start_th, end_th, dt,
1862 inputrec->opts.nFreeze, md->invmass, md->ptype,
1863 md->cFREEZE, md->cTC,
1864 x_rvec, xp_rvec, v_rvec, f_rvec,
1866 upd->sd->bd_rf.data(),
1867 step, inputrec->ld_seed, DOMAINDECOMP(cr) ? cr->dd->globalAtomIndices.data() : nullptr);
1872 gmx_bool bExtended = (inputrec->etc == etcNOSEHOOVER ||
1873 inputrec->epc == epcPARRINELLORAHMAN ||
1874 inputrec->epc == epcMTTK);
1876 /* assuming barostat coupled to group 0 */
1877 real alpha = 1.0 + DIM/static_cast<real>(inputrec->opts.nrdf[0]);
1882 do_update_vv_vel(start_th, end_th, dt,
1883 inputrec->opts.acc, inputrec->opts.nFreeze,
1884 md->invmass, md->ptype,
1885 md->cFREEZE, md->cACC,
1887 bExtended, state->veta, alpha);
1890 do_update_vv_pos(start_th, end_th, dt,
1891 inputrec->opts.nFreeze,
1892 md->ptype, md->cFREEZE,
1893 x_rvec, xp_rvec, v_rvec,
1894 bExtended, state->veta);
1900 gmx_fatal(FARGS, "Don't know how to update coordinates");
1903 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
1908 extern gmx_bool update_randomize_velocities(const t_inputrec *ir, int64_t step, const t_commrec *cr,
1909 const t_mdatoms *md,
1910 gmx::ArrayRef<gmx::RVec> v,
1911 const gmx_update_t *upd,
1912 const gmx::Constraints *constr)
1915 real rate = (ir->delta_t)/ir->opts.tau_t[0];
1917 if (ir->etc == etcANDERSEN && constr != nullptr)
1919 /* Currently, Andersen thermostat does not support constrained
1920 systems. Functionality exists in the andersen_tcoupl
1921 function in GROMACS 4.5.7 to allow this combination. That
1922 code could be ported to the current random-number
1923 generation approach, but has not yet been done because of
1924 lack of time and resources. */
1925 gmx_fatal(FARGS, "Normal Andersen is currently not supported with constraints, use massive Andersen instead");
1928 /* proceed with andersen if 1) it's fixed probability per
1929 particle andersen or 2) it's massive andersen and it's tau_t/dt */
1930 if ((ir->etc == etcANDERSEN) || do_per_step(step, roundToInt(1.0/rate)))
1932 andersen_tcoupl(ir, step, cr, md, v, rate,
1933 upd->sd->randomize_group,